Removed reference to nonexistent file in g_helix.
[gromacs/AngularHB.git] / include / pme.h
blob7ba5a1a366a15119970c8ceefa892be47e8a9c04
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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, 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.
39 #ifndef _pme_h
40 #define _pme_h
42 #include <stdio.h>
43 #include "visibility.h"
44 #include "typedefs.h"
45 #include "gmxcomplex.h"
46 #include "gmx_wallcycle.h"
48 #ifdef __cplusplus
49 extern "C" {
50 #endif
52 typedef real *splinevec[DIM];
54 enum { GMX_SUM_QGRID_FORWARD, GMX_SUM_QGRID_BACKWARD };
56 GMX_LIBMD_EXPORT
57 int gmx_pme_init(gmx_pme_t *pmedata,t_commrec *cr,
58 int nnodes_major,int nnodes_minor,
59 t_inputrec *ir,int homenr,
60 gmx_bool bFreeEnergy, gmx_bool bReproducible, int nthread);
61 /* Initialize the pme data structures resepectively.
62 * Return value 0 indicates all well, non zero is an error code.
65 GMX_LIBMD_EXPORT
66 int gmx_pme_reinit(gmx_pme_t * pmedata,
67 t_commrec * cr,
68 gmx_pme_t pme_src,
69 const t_inputrec * ir,
70 ivec grid_size);
71 /* As gmx_pme_init, but takes most settings, except the grid, from pme_src */
73 int gmx_pme_destroy(FILE *log,gmx_pme_t *pmedata);
74 /* Destroy the pme data structures resepectively.
75 * Return value 0 indicates all well, non zero is an error code.
78 #define GMX_PME_SPREAD_Q (1<<0)
79 #define GMX_PME_SOLVE (1<<1)
80 #define GMX_PME_CALC_F (1<<2)
81 #define GMX_PME_CALC_ENER_VIR (1<<3)
82 /* This forces the grid to be backtransformed even without GMX_PME_CALC_F */
83 #define GMX_PME_CALC_POT (1<<4)
84 #define GMX_PME_DO_ALL_F (GMX_PME_SPREAD_Q | GMX_PME_SOLVE | GMX_PME_CALC_F)
86 int gmx_pme_do(gmx_pme_t pme,
87 int start, int homenr,
88 rvec x[], rvec f[],
89 real chargeA[], real chargeB[],
90 matrix box, t_commrec *cr,
91 int maxshift_x, int maxshift_y,
92 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
93 matrix lrvir, real ewaldcoeff,
94 real *energy, real lambda,
95 real *dvdlambda, int flags);
96 /* Do a PME calculation for the long range electrostatics.
97 * flags, defined above, determine which parts of the calculation are performed.
98 * Return value 0 indicates all well, non zero is an error code.
101 GMX_LIBMD_EXPORT
102 int gmx_pmeonly(gmx_pme_t pme,
103 t_commrec *cr, t_nrnb *mynrnb,
104 gmx_wallcycle_t wcycle,
105 real ewaldcoeff, gmx_bool bGatherOnly,
106 t_inputrec *ir);
107 /* Called on the nodes that do PME exclusively (as slaves)
110 void gmx_pme_calc_energy(gmx_pme_t pme,int n,rvec *x,real *q,real *V);
111 /* Calculate the PME grid energy V for n charges with a potential
112 * in the pme struct determined before with a call to gmx_pme_do
113 * with at least GMX_PME_SPREAD_Q and GMX_PME_SOLVE specified.
114 * Note that the charges are not spread on the grid in the pme struct.
115 * Currently does not work in parallel or with free energy.
118 /* The following three routines are for PME/PP node splitting in pme_pp.c */
120 /* Abstract type for PME <-> PP communication */
121 typedef struct gmx_pme_pp *gmx_pme_pp_t;
123 gmx_pme_pp_t gmx_pme_pp_init(t_commrec *cr);
124 /* Initialize the PME-only side of the PME <-> PP communication */
126 void gmx_pme_send_q(t_commrec *cr,
127 gmx_bool bFreeEnergy, real *chargeA, real *chargeB,
128 int maxshift_x, int maxshift_y);
129 /* Send the charges and maxshift to out PME-only node. */
131 void gmx_pme_send_x(t_commrec *cr, matrix box, rvec *x,
132 gmx_bool bFreeEnergy, real lambda,
133 gmx_bool bEnerVir,
134 gmx_large_int_t step);
135 /* Send the coordinates to our PME-only node and request a PME calculation */
137 GMX_LIBMD_EXPORT
138 void gmx_pme_send_finish(t_commrec *cr);
139 /* Tell our PME-only node to finish */
141 GMX_LIBMD_EXPORT
142 void gmx_pme_send_switch(t_commrec *cr, ivec grid_size, real ewaldcoeff);
143 /* Tell our PME-only node to switch to a new grid size */
145 void gmx_pme_receive_f(t_commrec *cr,
146 rvec f[], matrix vir,
147 real *energy, real *dvdlambda,
148 float *pme_cycles);
149 /* PP nodes receive the long range forces from the PME nodes */
151 int gmx_pme_recv_q_x(gmx_pme_pp_t pme_pp,
152 real **chargeA, real **chargeB,
153 matrix box, rvec **x,rvec **f,
154 int *maxshift_x, int *maxshift_y,
155 gmx_bool *bFreeEnergy, real *lambda,
156 gmx_bool *bEnerVir,
157 gmx_large_int_t *step,
158 ivec grid_size, real *ewaldcoeff);
160 /* Receive charges and/or coordinates from the PP-only nodes.
161 * Returns the number of atoms, or -1 when the run is finished.
162 * In the special case of a PME grid size switch request, -2 is returned
163 * and grid_size and *ewaldcoeff are set, which are otherwise not set.
166 void gmx_pme_send_force_vir_ener(gmx_pme_pp_t pme_pp,
167 rvec *f, matrix vir,
168 real energy, real dvdlambda,
169 float cycles);
170 /* Send the PME mesh force, virial and energy to the PP-only nodes */
172 #ifdef __cplusplus
174 #endif
176 #endif