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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "visibility.h"
45 #include "gmxcomplex.h"
46 #include "gmx_wallcycle.h"
53 typedef real *splinevec[DIM];
56 GMX_SUM_QGRID_FORWARD, GMX_SUM_QGRID_BACKWARD
60 int gmx_pme_init(gmx_pme_t *pmedata, t_commrec *cr,
61 int nnodes_major, int nnodes_minor,
62 t_inputrec *ir, int homenr,
63 gmx_bool bFreeEnergy, gmx_bool bReproducible, int nthread);
64 /* Initialize the pme data structures resepectively.
65 * Return value 0 indicates all well, non zero is an error code.
69 int gmx_pme_reinit(gmx_pme_t * pmedata,
72 const t_inputrec * ir,
74 /* As gmx_pme_init, but takes most settings, except the grid, from pme_src */
76 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata);
77 /* Destroy the pme data structures resepectively.
78 * Return value 0 indicates all well, non zero is an error code.
81 #define GMX_PME_SPREAD_Q (1<<0)
82 #define GMX_PME_SOLVE (1<<1)
83 #define GMX_PME_CALC_F (1<<2)
84 #define GMX_PME_CALC_ENER_VIR (1<<3)
85 /* This forces the grid to be backtransformed even without GMX_PME_CALC_F */
86 #define GMX_PME_CALC_POT (1<<4)
87 #define GMX_PME_DO_ALL_F (GMX_PME_SPREAD_Q | GMX_PME_SOLVE | GMX_PME_CALC_F)
89 int gmx_pme_do(gmx_pme_t pme,
90 int start, int homenr,
92 real chargeA[], real chargeB[],
93 matrix box, t_commrec *cr,
94 int maxshift_x, int maxshift_y,
95 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
96 matrix lrvir, real ewaldcoeff,
97 real *energy, real lambda,
98 real *dvdlambda, int flags);
99 /* Do a PME calculation for the long range electrostatics.
100 * flags, defined above, determine which parts of the calculation are performed.
101 * Return value 0 indicates all well, non zero is an error code.
105 int gmx_pmeonly(gmx_pme_t pme,
106 t_commrec *cr, t_nrnb *mynrnb,
107 gmx_wallcycle_t wcycle,
108 gmx_runtime_t *runtime,
109 real ewaldcoeff, gmx_bool bGatherOnly,
111 /* Called on the nodes that do PME exclusively (as slaves)
114 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V);
115 /* Calculate the PME grid energy V for n charges with a potential
116 * in the pme struct determined before with a call to gmx_pme_do
117 * with at least GMX_PME_SPREAD_Q and GMX_PME_SOLVE specified.
118 * Note that the charges are not spread on the grid in the pme struct.
119 * Currently does not work in parallel or with free energy.
122 /* The following three routines are for PME/PP node splitting in pme_pp.c */
124 /* Abstract type for PME <-> PP communication */
125 typedef struct gmx_pme_pp *gmx_pme_pp_t;
128 void gmx_pme_check_restrictions(int pme_order,
129 int nkx, int nky, int nkz,
132 gmx_bool bUseThreads,
134 gmx_bool *bValidSettings);
135 /* Check restrictions on pme_order and the PME grid nkx,nky,nkz.
136 * With bFatal=TRUE, a fatal error is generated on violation,
137 * bValidSettings=NULL can be passed.
138 * With bFatal=FALSE, *bValidSettings reports the validity of the settings.
139 * bUseThreads tells if any MPI rank doing PME uses more than 1 threads.
140 * If at calling you bUseThreads is unknown, pass TRUE for conservative
144 gmx_pme_pp_t gmx_pme_pp_init(t_commrec *cr);
145 /* Initialize the PME-only side of the PME <-> PP communication */
147 void gmx_pme_send_q(t_commrec *cr,
148 gmx_bool bFreeEnergy, real *chargeA, real *chargeB,
149 int maxshift_x, int maxshift_y);
150 /* Send the charges and maxshift to out PME-only node. */
152 void gmx_pme_send_x(t_commrec *cr, matrix box, rvec *x,
153 gmx_bool bFreeEnergy, real lambda,
155 gmx_large_int_t step);
156 /* Send the coordinates to our PME-only node and request a PME calculation */
159 void gmx_pme_send_finish(t_commrec *cr);
160 /* Tell our PME-only node to finish */
163 void gmx_pme_send_switchgrid(t_commrec *cr, ivec grid_size, real ewaldcoeff);
164 /* Tell our PME-only node to switch to a new grid size */
167 void gmx_pme_send_resetcounters(t_commrec *cr, gmx_large_int_t step);
168 /* Tell our PME-only node to reset all cycle and flop counters */
170 void gmx_pme_receive_f(t_commrec *cr,
171 rvec f[], matrix vir,
172 real *energy, real *dvdlambda,
174 /* PP nodes receive the long range forces from the PME nodes */
176 /* Return values for gmx_pme_recv_q_x */
178 pmerecvqxX, /* calculate PME mesh interactions for new x */
179 pmerecvqxFINISH, /* the simulation should finish, we should quit */
180 pmerecvqxSWITCHGRID, /* change the PME grid size */
181 pmerecvqxRESETCOUNTERS /* reset the cycle and flop counters */
184 int gmx_pme_recv_q_x(gmx_pme_pp_t pme_pp,
186 real **chargeA, real **chargeB,
187 matrix box, rvec **x, rvec **f,
188 int *maxshift_x, int *maxshift_y,
189 gmx_bool *bFreeEnergy, real *lambda,
191 gmx_large_int_t *step,
192 ivec grid_size, real *ewaldcoeff);
194 /* With return value:
195 * pmerecvqxX: all parameters set, chargeA and chargeB can be NULL
196 * pmerecvqxFINISH: no parameters set
197 * pmerecvqxSWITCHGRID: only grid_size and *ewaldcoeff are set
198 * pmerecvqxRESETCOUNTERS: *step is set
201 void gmx_pme_send_force_vir_ener(gmx_pme_pp_t pme_pp,
203 real energy, real dvdlambda,
205 /* Send the PME mesh force, virial and energy to the PP-only nodes */