Merge release-5-0 into master
[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 "gmxpre.h"
36
37 #include "nbnxn_kernel_ref.h"
38
39 #include "config.h"
40
41 #include <assert.h>
42 #include <math.h>
43
44 #include "gromacs/legacyheaders/force.h"
45 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/nb_verlet.h"
49 #include "gromacs/mdlib/nbnxn_consts.h"
50 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_common.h"
51 #include "gromacs/pbcutil/ishift.h"
52 #include "gromacs/utility/smalloc.h"
53
54 /*! \brief Typedefs for declaring lookup tables of kernel functions.
55  */
56
57 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t     *nbl,
58                                   const nbnxn_atomdata_t     *nbat,
59                                   const interaction_const_t  *ic,
60                                   rvec                       *shift_vec,
61                                   real                       *f,
62                                   real                       *fshift);
63
64 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
65                                 const nbnxn_atomdata_t     *nbat,
66                                 const interaction_const_t  *ic,
67                                 rvec                       *shift_vec,
68                                 real                       *f,
69                                 real                       *fshift,
70                                 real                       *Vvdw,
71                                 real                       *Vc);
72
73 /* Analytical reaction-field kernels */
74 #define CALC_COUL_RF
75 #define LJ_CUT
76 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
77 #undef LJ_CUT
78 #define LJ_FORCE_SWITCH
79 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
80 #undef LJ_FORCE_SWITCH
81 #define LJ_POT_SWITCH
82 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
83 #undef LJ_POT_SWITCH
84 #define LJ_EWALD
85 #define LJ_CUT
86 #define LJ_EWALD_COMB_GEOM
87 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
88 #undef LJ_EWALD_COMB_GEOM
89 #define LJ_EWALD_COMB_LB
90 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
91 #undef LJ_EWALD_COMB_LB
92 #undef LJ_CUT
93 #undef LJ_EWALD
94 #undef CALC_COUL_RF
95
96
97 /* Tabulated exclusion interaction electrostatics kernels */
98 #define CALC_COUL_TAB
99 #define LJ_CUT
100 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
101 #undef LJ_CUT
102 #define LJ_FORCE_SWITCH
103 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
104 #undef LJ_FORCE_SWITCH
105 #define LJ_POT_SWITCH
106 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
107 #undef LJ_POT_SWITCH
108 #define LJ_EWALD
109 #define LJ_CUT
110 #define LJ_EWALD_COMB_GEOM
111 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
112 #undef LJ_EWALD_COMB_GEOM
113 #define LJ_EWALD_COMB_LB
114 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
115 #undef LJ_EWALD_COMB_LB
116 #undef LJ_CUT
117 #undef LJ_EWALD
118 /* Twin-range cut-off kernels */
119 #define VDW_CUTOFF_CHECK
120 #define LJ_CUT
121 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
122 #undef LJ_CUT
123 #define LJ_FORCE_SWITCH
124 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
125 #undef LJ_FORCE_SWITCH
126 #define LJ_POT_SWITCH
127 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
128 #undef LJ_POT_SWITCH
129 #define LJ_EWALD
130 #define LJ_CUT
131 #define LJ_EWALD_COMB_GEOM
132 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
133 #undef LJ_EWALD_COMB_GEOM
134 #define LJ_EWALD_COMB_LB
135 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref_includes.h"
136 #undef LJ_EWALD_COMB_LB
137 #undef LJ_CUT
138 #undef LJ_EWALD
139 #undef VDW_CUTOFF_CHECK
140 #undef CALC_COUL_TAB
141
142
143 enum {
144     coultRF, coultTAB, coultTAB_TWIN, coultNR
145 };
146
147 enum {
148     vdwtCUT, vdwtFSWITCH, vdwtPSWITCH, vdwtEWALDGEOM, vdwtEWALDLB, vdwtNR
149 };
150
151 p_nbk_func_noener p_nbk_c_noener[coultNR][vdwtNR] =
152 {
153     { 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           },
154     { 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        },
155     { 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 }
156 };
157
158 p_nbk_func_ener p_nbk_c_ener[coultNR][vdwtNR] =
159 {
160     { 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            },
161     { 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         },
162     { 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  }
163 };
164
165 p_nbk_func_ener p_nbk_c_energrp[coultNR][vdwtNR] =
166 {
167     { 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           },
168     { 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        },
169     { 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 }
170 };
171
172 void
173 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
174                  const nbnxn_atomdata_t     *nbat,
175                  const interaction_const_t  *ic,
176                  rvec                       *shift_vec,
177                  int                         force_flags,
178                  int                         clearF,
179                  real                       *fshift,
180                  real                       *Vc,
181                  real                       *Vvdw)
182 {
183     int                nnbl;
184     nbnxn_pairlist_t **nbl;
185     int                coult;
186     int                vdwt;
187     int                nb;
188     int                nthreads gmx_unused;
189
190     nnbl = nbl_list->nnbl;
191     nbl  = nbl_list->nbl;
192
193     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
194     {
195         coult = coultRF;
196     }
197     else
198     {
199         if (ic->rcoulomb == ic->rvdw)
200         {
201             coult = coultTAB;
202         }
203         else
204         {
205             coult = coultTAB_TWIN;
206         }
207     }
208
209     if (ic->vdwtype == evdwCUT)
210     {
211         switch (ic->vdw_modifier)
212         {
213             case eintmodPOTSHIFT:
214             case eintmodNONE:
215                 vdwt = vdwtCUT;
216                 break;
217             case eintmodFORCESWITCH:
218                 vdwt = vdwtFSWITCH;
219                 break;
220             case eintmodPOTSWITCH:
221                 vdwt = vdwtPSWITCH;
222                 break;
223             default:
224                 gmx_incons("Unsupported VdW modifier");
225                 break;
226         }
227     }
228     else if (ic->vdwtype == evdwPME)
229     {
230         if (ic->ljpme_comb_rule == ljcrGEOM)
231         {
232             assert(nbat->comb_rule == ljcrGEOM);
233             vdwt = vdwtEWALDGEOM;
234         }
235         else
236         {
237             assert(nbat->comb_rule == ljcrLB);
238             vdwt = vdwtEWALDLB;
239         }
240     }
241     else
242     {
243         gmx_incons("Unsupported vdwtype in nbnxn reference kernel");
244     }
245
246     nthreads = gmx_omp_nthreads_get(emntNonbonded);
247 #pragma omp parallel for schedule(static) num_threads(nthreads)
248     for (nb = 0; nb < nnbl; nb++)
249     {
250         nbnxn_atomdata_output_t *out;
251         real                    *fshift_p;
252
253         out = &nbat->out[nb];
254
255         if (clearF == enbvClearFYes)
256         {
257             clear_f(nbat, nb, out->f);
258         }
259
260         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
261         {
262             fshift_p = fshift;
263         }
264         else
265         {
266             fshift_p = out->fshift;
267
268             if (clearF == enbvClearFYes)
269             {
270                 clear_fshift(fshift_p);
271             }
272         }
273
274         if (!(force_flags & GMX_FORCE_ENERGY))
275         {
276             /* Don't calculate energies */
277             p_nbk_c_noener[coult][vdwt](nbl[nb], nbat,
278                                         ic,
279                                         shift_vec,
280                                         out->f,
281                                         fshift_p);
282         }
283         else if (out->nV == 1)
284         {
285             /* No energy groups */
286             out->Vvdw[0] = 0;
287             out->Vc[0]   = 0;
288
289             p_nbk_c_ener[coult][vdwt](nbl[nb], nbat,
290                                       ic,
291                                       shift_vec,
292                                       out->f,
293                                       fshift_p,
294                                       out->Vvdw,
295                                       out->Vc);
296         }
297         else
298         {
299             /* Calculate energy group contributions */
300             int i;
301
302             for (i = 0; i < out->nV; i++)
303             {
304                 out->Vvdw[i] = 0;
305             }
306             for (i = 0; i < out->nV; i++)
307             {
308                 out->Vc[i] = 0;
309             }
310
311             p_nbk_c_energrp[coult][vdwt](nbl[nb], nbat,
312                                          ic,
313                                          shift_vec,
314                                          out->f,
315                                          fshift_p,
316                                          out->Vvdw,
317                                          out->Vc);
318         }
319     }
320
321     if (force_flags & GMX_FORCE_ENERGY)
322     {
323         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
324     }
325 }