Sort all includes in src/gromacs
[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
189     nnbl = nbl_list->nnbl;
190     nbl  = nbl_list->nbl;
191
192     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
193     {
194         coult = coultRF;
195     }
196     else
197     {
198         if (ic->rcoulomb == ic->rvdw)
199         {
200             coult = coultTAB;
201         }
202         else
203         {
204             coult = coultTAB_TWIN;
205         }
206     }
207
208     if (ic->vdwtype == evdwCUT)
209     {
210         switch (ic->vdw_modifier)
211         {
212             case eintmodPOTSHIFT:
213             case eintmodNONE:
214                 vdwt = vdwtCUT;
215                 break;
216             case eintmodFORCESWITCH:
217                 vdwt = vdwtFSWITCH;
218                 break;
219             case eintmodPOTSWITCH:
220                 vdwt = vdwtPSWITCH;
221                 break;
222             default:
223                 gmx_incons("Unsupported VdW modifier");
224                 break;
225         }
226     }
227     else if (ic->vdwtype == evdwPME)
228     {
229         if (ic->ljpme_comb_rule == ljcrGEOM)
230         {
231             assert(nbat->comb_rule == ljcrGEOM);
232             vdwt = vdwtEWALDGEOM;
233         }
234         else
235         {
236             assert(nbat->comb_rule == ljcrLB);
237             vdwt = vdwtEWALDLB;
238         }
239     }
240     else
241     {
242         gmx_incons("Unsupported vdwtype in nbnxn reference kernel");
243     }
244
245 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
246     for (nb = 0; nb < nnbl; nb++)
247     {
248         nbnxn_atomdata_output_t *out;
249         real                    *fshift_p;
250
251         out = &nbat->out[nb];
252
253         if (clearF == enbvClearFYes)
254         {
255             clear_f(nbat, nb, out->f);
256         }
257
258         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
259         {
260             fshift_p = fshift;
261         }
262         else
263         {
264             fshift_p = out->fshift;
265
266             if (clearF == enbvClearFYes)
267             {
268                 clear_fshift(fshift_p);
269             }
270         }
271
272         if (!(force_flags & GMX_FORCE_ENERGY))
273         {
274             /* Don't calculate energies */
275             p_nbk_c_noener[coult][vdwt](nbl[nb], nbat,
276                                         ic,
277                                         shift_vec,
278                                         out->f,
279                                         fshift_p);
280         }
281         else if (out->nV == 1)
282         {
283             /* No energy groups */
284             out->Vvdw[0] = 0;
285             out->Vc[0]   = 0;
286
287             p_nbk_c_ener[coult][vdwt](nbl[nb], nbat,
288                                       ic,
289                                       shift_vec,
290                                       out->f,
291                                       fshift_p,
292                                       out->Vvdw,
293                                       out->Vc);
294         }
295         else
296         {
297             /* Calculate energy group contributions */
298             int i;
299
300             for (i = 0; i < out->nV; i++)
301             {
302                 out->Vvdw[i] = 0;
303             }
304             for (i = 0; i < out->nV; i++)
305             {
306                 out->Vc[i] = 0;
307             }
308
309             p_nbk_c_energrp[coult][vdwt](nbl[nb], nbat,
310                                          ic,
311                                          shift_vec,
312                                          out->f,
313                                          fshift_p,
314                                          out->Vvdw,
315                                          out->Vc);
316         }
317     }
318
319     if (force_flags & GMX_FORCE_ENERGY)
320     {
321         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
322     }
323 }