Mark advanced variables/options that are not relevant for the OpenMM build
[gromacs/rigid-bodies.git] / include / pme.h
blob75c76443e8d16ac92593975b6289a732464666ca
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _pme_h
37 #define _pme_h
39 #include <stdio.h>
40 #include "typedefs.h"
41 #include "gmxcomplex.h"
42 #include "gmx_wallcycle.h"
44 #ifdef __cplusplus
45 extern "C" {
46 #endif
48 typedef real *splinevec[DIM];
50 enum { GMX_SUM_QGRID_FORWARD, GMX_SUM_QGRID_BACKWARD };
52 int gmx_pme_init(gmx_pme_t *pmedata,t_commrec *cr,
53 int nnodes_major,int nnodes_minor,
54 t_inputrec *ir,int homenr,
55 gmx_bool bFreeEnergy, gmx_bool bReproducible);
57 int gmx_pme_destroy(FILE *log,gmx_pme_t *pmedata);
58 /* Initialize and destroy the pme data structures resepectively.
59 * Return value 0 indicates all well, non zero is an error code.
62 #define GMX_PME_SPREAD_Q (1<<0)
63 #define GMX_PME_SOLVE (1<<1)
64 #define GMX_PME_CALC_F (1<<2)
65 #define GMX_PME_CALC_ENER_VIR (1<<3)
66 #define GMX_PME_DO_ALL_F (GMX_PME_SPREAD_Q | GMX_PME_SOLVE | GMX_PME_CALC_F)
68 int gmx_pme_do(gmx_pme_t pme,
69 int start, int homenr,
70 rvec x[], rvec f[],
71 real chargeA[], real chargeB[],
72 matrix box, t_commrec *cr,
73 int maxshift_x, int maxshift_y,
74 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
75 matrix lrvir, real ewaldcoeff,
76 real *energy, real lambda,
77 real *dvdlambda, int flags);
78 /* Do a PME calculation for the long range electrostatics.
79 * flags, defined above, determine which parts of the calculation are performed.
80 * Return value 0 indicates all well, non zero is an error code.
83 int gmx_pmeonly(gmx_pme_t pme,
84 t_commrec *cr, t_nrnb *mynrnb,
85 gmx_wallcycle_t wcycle,
86 real ewaldcoeff, gmx_bool bGatherOnly,
87 t_inputrec *ir);
88 /* Called on the nodes that do PME exclusively (as slaves)
91 void gmx_pme_calc_energy(gmx_pme_t pme,int n,rvec *x,real *q,real *V);
92 /* Calculate the PME grid energy V for n charges with a potential
93 * in the pme struct determined before with a call to gmx_pme_do
94 * with at least GMX_PME_SPREAD_Q and GMX_PME_SOLVE specified.
95 * Note that the charges are not spread on the grid in the pme struct.
96 * Currently does not work in parallel or with free energy.
99 /* The following three routines are for PME/PP node splitting in pme_pp.c */
101 /* Abstract type for PME <-> PP communication */
102 typedef struct gmx_pme_pp *gmx_pme_pp_t;
104 gmx_pme_pp_t gmx_pme_pp_init(t_commrec *cr);
105 /* Initialize the PME-only side of the PME <-> PP communication */
107 void gmx_pme_send_q(t_commrec *cr,
108 gmx_bool bFreeEnergy, real *chargeA, real *chargeB,
109 int maxshift_x, int maxshift_y);
110 /* Send the charges and maxshift to out PME-only node. */
112 void gmx_pme_send_x(t_commrec *cr, matrix box, rvec *x,
113 gmx_bool bFreeEnergy, real lambda,
114 gmx_bool bEnerVir,
115 gmx_large_int_t step);
116 /* Send the coordinates to our PME-only node and request a PME calculation */
118 void gmx_pme_finish(t_commrec *cr);
119 /* Tell our PME-only node to finish */
121 void gmx_pme_receive_f(t_commrec *cr,
122 rvec f[], matrix vir,
123 real *energy, real *dvdlambda,
124 float *pme_cycles);
125 /* PP nodes receive the long range forces from the PME nodes */
127 int gmx_pme_recv_q_x(gmx_pme_pp_t pme_pp,
128 real **chargeA, real **chargeB,
129 matrix box, rvec **x,rvec **f,
130 int *maxshift_x,int *maxshift_y,
131 gmx_bool *bFreeEnergy,real *lambda,
132 gmx_bool *bEnerVir,
133 gmx_large_int_t *step);
134 /* Receive charges and/or coordinates from the PP-only nodes.
135 * Returns the number of atoms, or -1 when the run is finished.
138 void gmx_pme_send_force_vir_ener(gmx_pme_pp_t pme_pp,
139 rvec *f, matrix vir,
140 real energy, real dvdlambda,
141 float cycles);
142 /* Send the PME mesh force, virial and energy to the PP-only nodes */
144 #ifdef __cplusplus
146 #endif
148 #endif