Avoid using function calls in OpenMP directives
[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 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40 #include <assert.h>
41
42 #include "typedefs.h"
43 #include "vec.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "force.h"
46 #include "gmx_omp_nthreads.h"
47 #include "nbnxn_kernel_ref.h"
48 #include "../nbnxn_consts.h"
49 #include "nbnxn_kernel_common.h"
50
51 /*! \brief Typedefs for declaring lookup tables of kernel functions.
52  */
53
54 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t     *nbl,
55                                   const nbnxn_atomdata_t     *nbat,
56                                   const interaction_const_t  *ic,
57                                   rvec                       *shift_vec,
58                                   real                       *f,
59                                   real                       *fshift);
60
61 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
62                                 const nbnxn_atomdata_t     *nbat,
63                                 const interaction_const_t  *ic,
64                                 rvec                       *shift_vec,
65                                 real                       *f,
66                                 real                       *fshift,
67                                 real                       *Vvdw,
68                                 real                       *Vc);
69
70 /* Analytical reaction-field kernels */
71 #define CALC_COUL_RF
72 #define LJ_CUT
73 #include "nbnxn_kernel_ref_includes.h"
74 #undef LJ_CUT
75 #define LJ_FORCE_SWITCH
76 #include "nbnxn_kernel_ref_includes.h"
77 #undef LJ_FORCE_SWITCH
78 #define LJ_POT_SWITCH
79 #include "nbnxn_kernel_ref_includes.h"
80 #undef LJ_POT_SWITCH
81 #define LJ_EWALD
82 #define LJ_CUT
83 #define LJ_EWALD_COMB_GEOM
84 #include "nbnxn_kernel_ref_includes.h"
85 #undef LJ_EWALD_COMB_GEOM
86 #define LJ_EWALD_COMB_LB
87 #include "nbnxn_kernel_ref_includes.h"
88 #undef LJ_EWALD_COMB_LB
89 #undef LJ_CUT
90 #undef LJ_EWALD
91 #undef CALC_COUL_RF
92
93
94 /* Tabulated exclusion interaction electrostatics kernels */
95 #define CALC_COUL_TAB
96 #define LJ_CUT
97 #include "nbnxn_kernel_ref_includes.h"
98 #undef LJ_CUT
99 #define LJ_FORCE_SWITCH
100 #include "nbnxn_kernel_ref_includes.h"
101 #undef LJ_FORCE_SWITCH
102 #define LJ_POT_SWITCH
103 #include "nbnxn_kernel_ref_includes.h"
104 #undef LJ_POT_SWITCH
105 #define LJ_EWALD
106 #define LJ_CUT
107 #define LJ_EWALD_COMB_GEOM
108 #include "nbnxn_kernel_ref_includes.h"
109 #undef LJ_EWALD_COMB_GEOM
110 #define LJ_EWALD_COMB_LB
111 #include "nbnxn_kernel_ref_includes.h"
112 #undef LJ_EWALD_COMB_LB
113 #undef LJ_CUT
114 #undef LJ_EWALD
115 /* Twin-range cut-off kernels */
116 #define VDW_CUTOFF_CHECK
117 #define LJ_CUT
118 #include "nbnxn_kernel_ref_includes.h"
119 #undef LJ_CUT
120 #define LJ_FORCE_SWITCH
121 #include "nbnxn_kernel_ref_includes.h"
122 #undef LJ_FORCE_SWITCH
123 #define LJ_POT_SWITCH
124 #include "nbnxn_kernel_ref_includes.h"
125 #undef LJ_POT_SWITCH
126 #define LJ_EWALD
127 #define LJ_CUT
128 #define LJ_EWALD_COMB_GEOM
129 #include "nbnxn_kernel_ref_includes.h"
130 #undef LJ_EWALD_COMB_GEOM
131 #define LJ_EWALD_COMB_LB
132 #include "nbnxn_kernel_ref_includes.h"
133 #undef LJ_EWALD_COMB_LB
134 #undef LJ_CUT
135 #undef LJ_EWALD
136 #undef VDW_CUTOFF_CHECK
137 #undef CALC_COUL_TAB
138
139
140 enum {
141     coultRF, coultTAB, coultTAB_TWIN, coultNR
142 };
143
144 enum {
145     vdwtCUT, vdwtFSWITCH, vdwtPSWITCH, vdwtEWALDGEOM, vdwtEWALDLB, vdwtNR
146 };
147
148 p_nbk_func_noener p_nbk_c_noener[coultNR][vdwtNR] =
149 {
150     { 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           },
151     { 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        },
152     { 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 }
153 };
154
155 p_nbk_func_ener p_nbk_c_ener[coultNR][vdwtNR] =
156 {
157     { 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            },
158     { 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         },
159     { 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  }
160 };
161
162 p_nbk_func_ener p_nbk_c_energrp[coultNR][vdwtNR] =
163 {
164     { 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           },
165     { 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        },
166     { 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 }
167 };
168
169 void
170 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
171                  const nbnxn_atomdata_t     *nbat,
172                  const interaction_const_t  *ic,
173                  rvec                       *shift_vec,
174                  int                         force_flags,
175                  int                         clearF,
176                  real                       *fshift,
177                  real                       *Vc,
178                  real                       *Vvdw)
179 {
180     int                nnbl;
181     nbnxn_pairlist_t **nbl;
182     int                coult;
183     int                vdwt;
184     int                nb;
185     int                nthreads gmx_unused;
186
187     nnbl = nbl_list->nnbl;
188     nbl  = nbl_list->nbl;
189
190     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
191     {
192         coult = coultRF;
193     }
194     else
195     {
196         if (ic->rcoulomb == ic->rvdw)
197         {
198             coult = coultTAB;
199         }
200         else
201         {
202             coult = coultTAB_TWIN;
203         }
204     }
205
206     if (ic->vdwtype == evdwCUT)
207     {
208         switch (ic->vdw_modifier)
209         {
210             case eintmodPOTSHIFT:
211             case eintmodNONE:
212                 vdwt = vdwtCUT;
213                 break;
214             case eintmodFORCESWITCH:
215                 vdwt = vdwtFSWITCH;
216                 break;
217             case eintmodPOTSWITCH:
218                 vdwt = vdwtPSWITCH;
219                 break;
220             default:
221                 gmx_incons("Unsupported VdW modifier");
222                 break;
223         }
224     }
225     else if (ic->vdwtype == evdwPME)
226     {
227         if (ic->ljpme_comb_rule == ljcrGEOM)
228         {
229             assert(nbat->comb_rule == ljcrGEOM);
230             vdwt = vdwtEWALDGEOM;
231         }
232         else
233         {
234             assert(nbat->comb_rule == ljcrLB);
235             vdwt = vdwtEWALDLB;
236         }
237     }
238     else
239     {
240         gmx_incons("Unsupported vdwtype in nbnxn reference kernel");
241     }
242
243     nthreads = gmx_omp_nthreads_get(emntNonbonded);
244 #pragma omp parallel for schedule(static) num_threads(nthreads)
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 }