Remove all unnecessary HAVE_CONFIG_H
[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
52 /*! \brief Typedefs for declaring lookup tables of kernel functions.
53  */
54
55 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t     *nbl,
56                                   const nbnxn_atomdata_t     *nbat,
57                                   const interaction_const_t  *ic,
58                                   rvec                       *shift_vec,
59                                   real                       *f,
60                                   real                       *fshift);
61
62 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
63                                 const nbnxn_atomdata_t     *nbat,
64                                 const interaction_const_t  *ic,
65                                 rvec                       *shift_vec,
66                                 real                       *f,
67                                 real                       *fshift,
68                                 real                       *Vvdw,
69                                 real                       *Vc);
70
71 /* Analytical reaction-field kernels */
72 #define CALC_COUL_RF
73 #define LJ_CUT
74 #include "nbnxn_kernel_ref_includes.h"
75 #undef LJ_CUT
76 #define LJ_FORCE_SWITCH
77 #include "nbnxn_kernel_ref_includes.h"
78 #undef LJ_FORCE_SWITCH
79 #define LJ_POT_SWITCH
80 #include "nbnxn_kernel_ref_includes.h"
81 #undef LJ_POT_SWITCH
82 #define LJ_EWALD
83 #define LJ_CUT
84 #define LJ_EWALD_COMB_GEOM
85 #include "nbnxn_kernel_ref_includes.h"
86 #undef LJ_EWALD_COMB_GEOM
87 #define LJ_EWALD_COMB_LB
88 #include "nbnxn_kernel_ref_includes.h"
89 #undef LJ_EWALD_COMB_LB
90 #undef LJ_CUT
91 #undef LJ_EWALD
92 #undef CALC_COUL_RF
93
94
95 /* Tabulated exclusion interaction electrostatics kernels */
96 #define CALC_COUL_TAB
97 #define LJ_CUT
98 #include "nbnxn_kernel_ref_includes.h"
99 #undef LJ_CUT
100 #define LJ_FORCE_SWITCH
101 #include "nbnxn_kernel_ref_includes.h"
102 #undef LJ_FORCE_SWITCH
103 #define LJ_POT_SWITCH
104 #include "nbnxn_kernel_ref_includes.h"
105 #undef LJ_POT_SWITCH
106 #define LJ_EWALD
107 #define LJ_CUT
108 #define LJ_EWALD_COMB_GEOM
109 #include "nbnxn_kernel_ref_includes.h"
110 #undef LJ_EWALD_COMB_GEOM
111 #define LJ_EWALD_COMB_LB
112 #include "nbnxn_kernel_ref_includes.h"
113 #undef LJ_EWALD_COMB_LB
114 #undef LJ_CUT
115 #undef LJ_EWALD
116 /* Twin-range cut-off kernels */
117 #define VDW_CUTOFF_CHECK
118 #define LJ_CUT
119 #include "nbnxn_kernel_ref_includes.h"
120 #undef LJ_CUT
121 #define LJ_FORCE_SWITCH
122 #include "nbnxn_kernel_ref_includes.h"
123 #undef LJ_FORCE_SWITCH
124 #define LJ_POT_SWITCH
125 #include "nbnxn_kernel_ref_includes.h"
126 #undef LJ_POT_SWITCH
127 #define LJ_EWALD
128 #define LJ_CUT
129 #define LJ_EWALD_COMB_GEOM
130 #include "nbnxn_kernel_ref_includes.h"
131 #undef LJ_EWALD_COMB_GEOM
132 #define LJ_EWALD_COMB_LB
133 #include "nbnxn_kernel_ref_includes.h"
134 #undef LJ_EWALD_COMB_LB
135 #undef LJ_CUT
136 #undef LJ_EWALD
137 #undef VDW_CUTOFF_CHECK
138 #undef CALC_COUL_TAB
139
140
141 enum {
142     coultRF, coultTAB, coultTAB_TWIN, coultNR
143 };
144
145 enum {
146     vdwtCUT, vdwtFSWITCH, vdwtPSWITCH, vdwtEWALDGEOM, vdwtEWALDLB, vdwtNR
147 };
148
149 p_nbk_func_noener p_nbk_c_noener[coultNR][vdwtNR] =
150 {
151     { 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           },
152     { 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        },
153     { 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 }
154 };
155
156 p_nbk_func_ener p_nbk_c_ener[coultNR][vdwtNR] =
157 {
158     { 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            },
159     { 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         },
160     { 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  }
161 };
162
163 p_nbk_func_ener p_nbk_c_energrp[coultNR][vdwtNR] =
164 {
165     { 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           },
166     { 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        },
167     { 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 }
168 };
169
170 void
171 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
172                  const nbnxn_atomdata_t     *nbat,
173                  const interaction_const_t  *ic,
174                  rvec                       *shift_vec,
175                  int                         force_flags,
176                  int                         clearF,
177                  real                       *fshift,
178                  real                       *Vc,
179                  real                       *Vvdw)
180 {
181     int                nnbl;
182     nbnxn_pairlist_t **nbl;
183     int                coult;
184     int                vdwt;
185     int                nb;
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 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
244     for (nb = 0; nb < nnbl; nb++)
245     {
246         nbnxn_atomdata_output_t *out;
247         real                    *fshift_p;
248
249         out = &nbat->out[nb];
250
251         if (clearF == enbvClearFYes)
252         {
253             clear_f(nbat, nb, out->f);
254         }
255
256         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
257         {
258             fshift_p = fshift;
259         }
260         else
261         {
262             fshift_p = out->fshift;
263
264             if (clearF == enbvClearFYes)
265             {
266                 clear_fshift(fshift_p);
267             }
268         }
269
270         if (!(force_flags & GMX_FORCE_ENERGY))
271         {
272             /* Don't calculate energies */
273             p_nbk_c_noener[coult][vdwt](nbl[nb], nbat,
274                                         ic,
275                                         shift_vec,
276                                         out->f,
277                                         fshift_p);
278         }
279         else if (out->nV == 1)
280         {
281             /* No energy groups */
282             out->Vvdw[0] = 0;
283             out->Vc[0]   = 0;
284
285             p_nbk_c_ener[coult][vdwt](nbl[nb], nbat,
286                                       ic,
287                                       shift_vec,
288                                       out->f,
289                                       fshift_p,
290                                       out->Vvdw,
291                                       out->Vc);
292         }
293         else
294         {
295             /* Calculate energy group contributions */
296             int i;
297
298             for (i = 0; i < out->nV; i++)
299             {
300                 out->Vvdw[i] = 0;
301             }
302             for (i = 0; i < out->nV; i++)
303             {
304                 out->Vc[i] = 0;
305             }
306
307             p_nbk_c_energrp[coult][vdwt](nbl[nb], nbat,
308                                          ic,
309                                          shift_vec,
310                                          out->f,
311                                          fshift_p,
312                                          out->Vvdw,
313                                          out->Vc);
314         }
315     }
316
317     if (force_flags & GMX_FORCE_ENERGY)
318     {
319         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
320     }
321 }