Move some verlet headers to mdlib
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_ref.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "config.h"
36
37 #include <math.h>
38 #include <assert.h>
39
40 #include "typedefs.h"
41 #include "force.h"
42 #include "gmx_omp_nthreads.h"
43
44 #include "gromacs/math/vec.h"
45 #include "gromacs/pbcutil/ishift.h"
46 #include "gromacs/utility/smalloc.h"
47
48 #include "nbnxn_kernel_ref.h"
49 #include "../nbnxn_consts.h"
50 #include "nbnxn_kernel_common.h"
51 #include "gromacs/mdlib/nb_verlet.h"
52
53 /*! \brief Typedefs for declaring lookup tables of kernel functions.
54  */
55
56 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t     *nbl,
57                                   const nbnxn_atomdata_t     *nbat,
58                                   const interaction_const_t  *ic,
59                                   rvec                       *shift_vec,
60                                   real                       *f,
61                                   real                       *fshift);
62
63 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
64                                 const nbnxn_atomdata_t     *nbat,
65                                 const interaction_const_t  *ic,
66                                 rvec                       *shift_vec,
67                                 real                       *f,
68                                 real                       *fshift,
69                                 real                       *Vvdw,
70                                 real                       *Vc);
71
72 /* Analytical reaction-field kernels */
73 #define CALC_COUL_RF
74 #define LJ_CUT
75 #include "nbnxn_kernel_ref_includes.h"
76 #undef LJ_CUT
77 #define LJ_FORCE_SWITCH
78 #include "nbnxn_kernel_ref_includes.h"
79 #undef LJ_FORCE_SWITCH
80 #define LJ_POT_SWITCH
81 #include "nbnxn_kernel_ref_includes.h"
82 #undef LJ_POT_SWITCH
83 #define LJ_EWALD
84 #define LJ_CUT
85 #define LJ_EWALD_COMB_GEOM
86 #include "nbnxn_kernel_ref_includes.h"
87 #undef LJ_EWALD_COMB_GEOM
88 #define LJ_EWALD_COMB_LB
89 #include "nbnxn_kernel_ref_includes.h"
90 #undef LJ_EWALD_COMB_LB
91 #undef LJ_CUT
92 #undef LJ_EWALD
93 #undef CALC_COUL_RF
94
95
96 /* Tabulated exclusion interaction electrostatics kernels */
97 #define CALC_COUL_TAB
98 #define LJ_CUT
99 #include "nbnxn_kernel_ref_includes.h"
100 #undef LJ_CUT
101 #define LJ_FORCE_SWITCH
102 #include "nbnxn_kernel_ref_includes.h"
103 #undef LJ_FORCE_SWITCH
104 #define LJ_POT_SWITCH
105 #include "nbnxn_kernel_ref_includes.h"
106 #undef LJ_POT_SWITCH
107 #define LJ_EWALD
108 #define LJ_CUT
109 #define LJ_EWALD_COMB_GEOM
110 #include "nbnxn_kernel_ref_includes.h"
111 #undef LJ_EWALD_COMB_GEOM
112 #define LJ_EWALD_COMB_LB
113 #include "nbnxn_kernel_ref_includes.h"
114 #undef LJ_EWALD_COMB_LB
115 #undef LJ_CUT
116 #undef LJ_EWALD
117 /* Twin-range cut-off kernels */
118 #define VDW_CUTOFF_CHECK
119 #define LJ_CUT
120 #include "nbnxn_kernel_ref_includes.h"
121 #undef LJ_CUT
122 #define LJ_FORCE_SWITCH
123 #include "nbnxn_kernel_ref_includes.h"
124 #undef LJ_FORCE_SWITCH
125 #define LJ_POT_SWITCH
126 #include "nbnxn_kernel_ref_includes.h"
127 #undef LJ_POT_SWITCH
128 #define LJ_EWALD
129 #define LJ_CUT
130 #define LJ_EWALD_COMB_GEOM
131 #include "nbnxn_kernel_ref_includes.h"
132 #undef LJ_EWALD_COMB_GEOM
133 #define LJ_EWALD_COMB_LB
134 #include "nbnxn_kernel_ref_includes.h"
135 #undef LJ_EWALD_COMB_LB
136 #undef LJ_CUT
137 #undef LJ_EWALD
138 #undef VDW_CUTOFF_CHECK
139 #undef CALC_COUL_TAB
140
141
142 enum {
143     coultRF, coultTAB, coultTAB_TWIN, coultNR
144 };
145
146 enum {
147     vdwtCUT, vdwtFSWITCH, vdwtPSWITCH, vdwtEWALDGEOM, vdwtEWALDLB, vdwtNR
148 };
149
150 p_nbk_func_noener p_nbk_c_noener[coultNR][vdwtNR] =
151 {
152     { nbnxn_kernel_ElecRF_VdwLJ_F_ref,           nbnxn_kernel_ElecRF_VdwLJFsw_F_ref,           nbnxn_kernel_ElecRF_VdwLJPsw_F_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_ref           },
153     { nbnxn_kernel_ElecQSTab_VdwLJ_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJFsw_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJPsw_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombLB_F_ref        },
154     { nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJFsw_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJPsw_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombLB_F_ref }
155 };
156
157 p_nbk_func_ener p_nbk_c_ener[coultNR][vdwtNR] =
158 {
159     { nbnxn_kernel_ElecRF_VdwLJ_VF_ref,           nbnxn_kernel_ElecRF_VdwLJFsw_VF_ref,           nbnxn_kernel_ElecRF_VdwLJPsw_VF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_ref            },
160     { nbnxn_kernel_ElecQSTab_VdwLJ_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJFsw_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJPsw_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombLB_VF_ref         },
161     { nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJFsw_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJPsw_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombLB_VF_ref  }
162 };
163
164 p_nbk_func_ener p_nbk_c_energrp[coultNR][vdwtNR] =
165 {
166     { nbnxn_kernel_ElecRF_VdwLJ_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJFsw_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJPsw_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombLB_VgrpF_ref           },
167     { nbnxn_kernel_ElecQSTab_VdwLJ_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJFsw_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJPsw_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombLB_VgrpF_ref        },
168     { nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJFsw_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJPsw_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombLB_VgrpF_ref }
169 };
170
171 void
172 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
173                  const nbnxn_atomdata_t     *nbat,
174                  const interaction_const_t  *ic,
175                  rvec                       *shift_vec,
176                  int                         force_flags,
177                  int                         clearF,
178                  real                       *fshift,
179                  real                       *Vc,
180                  real                       *Vvdw)
181 {
182     int                nnbl;
183     nbnxn_pairlist_t **nbl;
184     int                coult;
185     int                vdwt;
186     int                nb;
187
188     nnbl = nbl_list->nnbl;
189     nbl  = nbl_list->nbl;
190
191     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
192     {
193         coult = coultRF;
194     }
195     else
196     {
197         if (ic->rcoulomb == ic->rvdw)
198         {
199             coult = coultTAB;
200         }
201         else
202         {
203             coult = coultTAB_TWIN;
204         }
205     }
206
207     if (ic->vdwtype == evdwCUT)
208     {
209         switch (ic->vdw_modifier)
210         {
211             case eintmodPOTSHIFT:
212             case eintmodNONE:
213                 vdwt = vdwtCUT;
214                 break;
215             case eintmodFORCESWITCH:
216                 vdwt = vdwtFSWITCH;
217                 break;
218             case eintmodPOTSWITCH:
219                 vdwt = vdwtPSWITCH;
220                 break;
221             default:
222                 gmx_incons("Unsupported VdW modifier");
223                 break;
224         }
225     }
226     else if (ic->vdwtype == evdwPME)
227     {
228         if (ic->ljpme_comb_rule == ljcrGEOM)
229         {
230             assert(nbat->comb_rule == ljcrGEOM);
231             vdwt = vdwtEWALDGEOM;
232         }
233         else
234         {
235             assert(nbat->comb_rule == ljcrLB);
236             vdwt = vdwtEWALDLB;
237         }
238     }
239     else
240     {
241         gmx_incons("Unsupported vdwtype in nbnxn reference kernel");
242     }
243
244 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
245     for (nb = 0; nb < nnbl; nb++)
246     {
247         nbnxn_atomdata_output_t *out;
248         real                    *fshift_p;
249
250         out = &nbat->out[nb];
251
252         if (clearF == enbvClearFYes)
253         {
254             clear_f(nbat, nb, out->f);
255         }
256
257         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
258         {
259             fshift_p = fshift;
260         }
261         else
262         {
263             fshift_p = out->fshift;
264
265             if (clearF == enbvClearFYes)
266             {
267                 clear_fshift(fshift_p);
268             }
269         }
270
271         if (!(force_flags & GMX_FORCE_ENERGY))
272         {
273             /* Don't calculate energies */
274             p_nbk_c_noener[coult][vdwt](nbl[nb], nbat,
275                                         ic,
276                                         shift_vec,
277                                         out->f,
278                                         fshift_p);
279         }
280         else if (out->nV == 1)
281         {
282             /* No energy groups */
283             out->Vvdw[0] = 0;
284             out->Vc[0]   = 0;
285
286             p_nbk_c_ener[coult][vdwt](nbl[nb], nbat,
287                                       ic,
288                                       shift_vec,
289                                       out->f,
290                                       fshift_p,
291                                       out->Vvdw,
292                                       out->Vc);
293         }
294         else
295         {
296             /* Calculate energy group contributions */
297             int i;
298
299             for (i = 0; i < out->nV; i++)
300             {
301                 out->Vvdw[i] = 0;
302             }
303             for (i = 0; i < out->nV; i++)
304             {
305                 out->Vc[i] = 0;
306             }
307
308             p_nbk_c_energrp[coult][vdwt](nbl[nb], nbat,
309                                          ic,
310                                          shift_vec,
311                                          out->f,
312                                          fshift_p,
313                                          out->Vvdw,
314                                          out->Vc);
315         }
316     }
317
318     if (force_flags & GMX_FORCE_ENERGY)
319     {
320         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
321     }
322 }