Move gmxcomplex.h to gromacs/math/
[alexxy/gromacs.git] / src / gromacs / legacyheaders / pme.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37
38 #ifndef _pme_h
39 #define _pme_h
40
41 #include <stdio.h>
42 #include "typedefs.h"
43 #include "sim_util.h"
44
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
48
49 typedef real *splinevec[DIM];
50
51 enum {
52     GMX_SUM_QGRID_FORWARD, GMX_SUM_QGRID_BACKWARD
53 };
54
55 int gmx_pme_init(gmx_pme_t *pmedata, t_commrec *cr,
56                  int nnodes_major, int nnodes_minor,
57                  t_inputrec *ir, int homenr,
58                  gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
59                  gmx_bool bReproducible, int nthread);
60 /* Initialize the pme data structures resepectively.
61  * Return value 0 indicates all well, non zero is an error code.
62  */
63
64 int gmx_pme_reinit(gmx_pme_t *         pmedata,
65                    t_commrec *         cr,
66                    gmx_pme_t           pme_src,
67                    const t_inputrec *  ir,
68                    ivec                grid_size);
69 /* As gmx_pme_init, but takes most settings, except the grid, from pme_src */
70
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.
74  */
75
76 #define GMX_PME_SPREAD_Q      (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)
82
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
85  */
86 #define GMX_PME_DO_COULOMB    (1<<13)
87 #define GMX_PME_DO_LJ         (1<<14)
88 #define GMX_PME_LJ_LB         (1<<15)
89
90 #define GMX_PME_DO_ALL_F  (GMX_PME_SPREAD_Q | GMX_PME_SOLVE | GMX_PME_CALC_F)
91
92 int gmx_pme_do(gmx_pme_t pme,
93                int start,       int homenr,
94                rvec x[],        rvec f[],
95                real chargeA[],  real chargeB[],
96                real c6A[],      real c6B[],
97                real sigmaA[],   real sigmaB[],
98                matrix box,      t_commrec *cr,
99                int  maxshift_x, int maxshift_y,
100                t_nrnb *nrnb,    gmx_wallcycle_t wcycle,
101                matrix vir_q,    real ewaldcoeff_q,
102                matrix vir_lj,   real ewaldcoeff_lj,
103                real *energy_q,  real *energy_lj,
104                real lambda_q, real lambda_lj,
105                real *dvdlambda_q, real *dvdlambda_lj,
106                int flags);
107 /* Do a PME calculation for the long range electrostatics and/or LJ.
108  * flags, defined above, determine which parts of the calculation are performed.
109  * Return value 0 indicates all well, non zero is an error code.
110  */
111
112 int gmx_pmeonly(gmx_pme_t pme,
113                 t_commrec *cr,     t_nrnb *mynrnb,
114                 gmx_wallcycle_t wcycle,
115                 gmx_walltime_accounting_t walltime_accounting,
116                 real ewaldcoeff_q, real ewaldcoeff_lj,
117                 t_inputrec *ir);
118 /* Called on the nodes that do PME exclusively (as slaves)
119  */
120
121 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V);
122 /* Calculate the PME grid energy V for n charges with a potential
123  * in the pme struct determined before with a call to gmx_pme_do
124  * with at least GMX_PME_SPREAD_Q and GMX_PME_SOLVE specified.
125  * Note that the charges are not spread on the grid in the pme struct.
126  * Currently does not work in parallel or with free energy.
127  */
128
129 /* The following three routines are for PME/PP node splitting in pme_pp.c */
130
131 /* Abstract type for PME <-> PP communication */
132 typedef struct gmx_pme_pp *gmx_pme_pp_t;
133
134 void gmx_pme_check_restrictions(int pme_order,
135                                 int nkx, int nky, int nkz,
136                                 int nnodes_major,
137                                 int nnodes_minor,
138                                 gmx_bool bUseThreads,
139                                 gmx_bool bFatal,
140                                 gmx_bool *bValidSettings);
141 /* Check restrictions on pme_order and the PME grid nkx,nky,nkz.
142  * With bFatal=TRUE, a fatal error is generated on violation,
143  * bValidSettings=NULL can be passed.
144  * With bFatal=FALSE, *bValidSettings reports the validity of the settings.
145  * bUseThreads tells if any MPI rank doing PME uses more than 1 threads.
146  * If at calling you bUseThreads is unknown, pass TRUE for conservative
147  * checking.
148  */
149
150 gmx_pme_pp_t gmx_pme_pp_init(t_commrec *cr);
151 /* Initialize the PME-only side of the PME <-> PP communication */
152
153 void gmx_pme_send_parameters(t_commrec *cr,
154                              gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
155                              real *chargeA, real *chargeB,
156                              real *c6A, real *c6B, real *sigmaA, real *sigmaB,
157                              int maxshift_x, int maxshift_y);
158 /* Send the charges and maxshift to out PME-only node. */
159
160 void gmx_pme_send_coordinates(t_commrec *cr, matrix box, rvec *x,
161                               gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
162                               real lambda_q, real lambda_lj,
163                               gmx_bool bEnerVir, int pme_flags,
164                               gmx_int64_t step);
165 /* Send the coordinates to our PME-only node and request a PME calculation */
166
167 void gmx_pme_send_finish(t_commrec *cr);
168 /* Tell our PME-only node to finish */
169
170 void gmx_pme_send_switchgrid(t_commrec *cr, ivec grid_size, real ewaldcoeff_q, real ewaldcoeff_lj);
171 /* Tell our PME-only node to switch to a new grid size */
172
173 void gmx_pme_send_resetcounters(t_commrec *cr, gmx_int64_t step);
174 /* Tell our PME-only node to reset all cycle and flop counters */
175
176 void gmx_pme_receive_f(t_commrec *cr,
177                        rvec f[], matrix vir_q, real *energy_q,
178                        matrix vir_lj, real *energy_lj,
179                        real *dvdlambda_q, real *dvdlambda_lj,
180                        float *pme_cycles);
181 /* PP nodes receive the long range forces from the PME nodes */
182
183 /* Return values for gmx_pme_recv_q_x */
184 enum {
185     pmerecvqxX,            /* calculate PME mesh interactions for new x    */
186     pmerecvqxFINISH,       /* the simulation should finish, we should quit */
187     pmerecvqxSWITCHGRID,   /* change the PME grid size                     */
188     pmerecvqxRESETCOUNTERS /* reset the cycle and flop counters            */
189 };
190
191 int gmx_pme_recv_params_coords(gmx_pme_pp_t pme_pp,
192                                int *natoms,
193                                real **chargeA, real **chargeB,
194                                real **c6A, real **c6B,
195                                real **sigmaA, real **sigmaB,
196                                matrix box, rvec **x, rvec **f,
197                                int *maxshift_x, int *maxshift_y,
198                                gmx_bool *bFreeEnergy_q, gmx_bool *bFreeEnergy_lj,
199                                real *lambda_q, real *lambda_lj,
200                                gmx_bool *bEnerVir, int *pme_flags,
201                                gmx_int64_t *step,
202                                ivec grid_size, real *ewaldcoeff_q, real *ewaldcoeff_lj);
203 ;
204 /* With return value:
205  * pmerecvqxX:             all parameters set, chargeA and chargeB can be NULL
206  * pmerecvqxFINISH:        no parameters set
207  * pmerecvqxSWITCHGRID:    only grid_size and *ewaldcoeff are set
208  * pmerecvqxRESETCOUNTERS: *step is set
209  */
210
211 void gmx_pme_send_force_vir_ener(gmx_pme_pp_t pme_pp,
212                                  rvec *f, matrix vir_q, real energy_q,
213                                  matrix vir_lj, real energy_lj,
214                                  real dvdlambda_q, real dvdlambda_lj,
215                                  float cycles);
216 /* Send the PME mesh force, virial and energy to the PP-only nodes */
217
218 #ifdef __cplusplus
219 }
220 #endif
221
222 #endif