Implemented nbnxn LJ switch functions
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_file_generator / nbnxn_kernel_simd_template.c.pre
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 "typedefs.h"
40
41 #ifdef {0}
42
43 {1}
44 #include "gromacs/simd/macros.h"
45 #include "gromacs/simd/vector_operations.h"
46 {2}
47 #define GMX_SIMD_J_UNROLL_SIZE {3}
48 #include "{4}"
49 #include "../nbnxn_kernel_common.h"
50 #include "gmx_omp_nthreads.h"
51 #include "types/force_flags.h"
52 #include "gmx_fatal.h"
53
54 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
55  */
56 enum {{
57     coultRF, coultTAB, coultTAB_TWIN, coultEWALD, coultEWALD_TWIN, coultNR
58 }};
59
60 /* Declare and define the kernel function pointer lookup tables.
61  * The minor index of the array goes over both the LJ combination rules,
62  * which is only supported by plain cut-off, and the LJ switch functions.
63  */
64 static p_nbk_func_noener p_nbk_noener[coultNR][ljcrNR+2] =
65 {7}
66 static p_nbk_func_ener p_nbk_ener[coultNR][ljcrNR+2] =
67 {8}
68 static p_nbk_func_ener p_nbk_energrp[coultNR][ljcrNR+2] =
69 {9}
70
71 static void
72 reduce_group_energies(int ng, int ng_2log,
73                       const real *VSvdw, const real *VSc,
74                       real *Vvdw, real *Vc)
75 {{
76     const int unrollj      = GMX_SIMD_REAL_WIDTH/GMX_SIMD_J_UNROLL_SIZE;
77     const int unrollj_half = unrollj/2;
78     int       ng_p2, i, j, j0, j1, c, s;
79
80     ng_p2 = (1<<ng_2log);
81
82     /* The size of the x86 SIMD energy group buffer array is:
83      * ng*ng*ng_p2*unrollj_half*simd_width
84      */
85     for (i = 0; i < ng; i++)
86     {{
87         for (j = 0; j < ng; j++)
88         {{
89             Vvdw[i*ng+j] = 0;
90             Vc[i*ng+j]   = 0;
91         }}
92
93         for (j1 = 0; j1 < ng; j1++)
94         {{
95             for (j0 = 0; j0 < ng; j0++)
96             {{
97                 c = ((i*ng + j1)*ng_p2 + j0)*unrollj_half*unrollj;
98                 for (s = 0; s < unrollj_half; s++)
99                 {{
100                     Vvdw[i*ng+j0] += VSvdw[c+0];
101                     Vvdw[i*ng+j1] += VSvdw[c+1];
102                     Vc  [i*ng+j0] += VSc  [c+0];
103                     Vc  [i*ng+j1] += VSc  [c+1];
104                     c             += unrollj + 2;
105                 }}
106             }}
107         }}
108     }}
109 }}
110
111 #else /* {0} */
112
113 #include "gmx_fatal.h"
114
115 #endif /* {0} */
116
117 void
118 {5}(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
119 {6}const nbnxn_atomdata_t    gmx_unused *nbat,
120 {6}const interaction_const_t gmx_unused *ic,
121 {6}int                       gmx_unused  ewald_excl,
122 {6}rvec                      gmx_unused *shift_vec,
123 {6}int                       gmx_unused  force_flags,
124 {6}int                       gmx_unused  clearF,
125 {6}real                      gmx_unused *fshift,
126 {6}real                      gmx_unused *Vc,
127 {6}real                      gmx_unused *Vvdw)
128 #ifdef {0}
129 {{
130     int                nnbl;
131     nbnxn_pairlist_t **nbl;
132     int                coult, ljtreatment = 0;
133     int                nb;
134
135     nnbl = nbl_list->nnbl;
136     nbl  = nbl_list->nbl;
137
138     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
139     {{
140         coult = coultRF;
141     }}
142     else
143     {{
144         if (ewald_excl == ewaldexclTable)
145         {{
146             if (ic->rcoulomb == ic->rvdw)
147             {{
148                 coult = coultTAB;
149             }}
150             else
151             {{
152                 coult = coultTAB_TWIN;
153             }}
154         }}
155         else
156         {{
157             if (ic->rcoulomb == ic->rvdw)
158             {{
159                 coult = coultEWALD;
160             }}
161             else
162             {{
163                 coult = coultEWALD_TWIN;
164             }}
165         }}
166     }}
167
168     switch (ic->vdw_modifier)
169     {{
170         case eintmodNONE:
171         case eintmodPOTSHIFT:
172             ljtreatment = nbat->comb_rule;
173             break;
174         /* Switch functions follow after cut-off combination rule kernels */
175         case eintmodFORCESWITCH:
176             ljtreatment = ljcrNR;
177             break;
178         case eintmodPOTSWITCH:
179             ljtreatment = ljcrNR + 1;
180             break;
181         default:
182             gmx_incons("Unsupported VdW interaction modifier");
183     }}
184
185 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
186     for (nb = 0; nb < nnbl; nb++)
187     {{
188         nbnxn_atomdata_output_t *out;
189         real                    *fshift_p;
190
191         out = &nbat->out[nb];
192
193         if (clearF == enbvClearFYes)
194         {{
195             clear_f(nbat, nb, out->f);
196         }}
197
198         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
199         {{
200             fshift_p = fshift;
201         }}
202         else
203         {{
204             fshift_p = out->fshift;
205
206             if (clearF == enbvClearFYes)
207             {{
208                 clear_fshift(fshift_p);
209             }}
210         }}
211
212         if (!(force_flags & GMX_FORCE_ENERGY))
213         {{
214             /* Don't calculate energies */
215             p_nbk_noener[coult][ljtreatment](nbl[nb], nbat,
216                                              ic,
217                                              shift_vec,
218                                              out->f,
219                                              fshift_p);
220         }}
221         else if (out->nV == 1)
222         {{
223             /* No energy groups */
224             out->Vvdw[0] = 0;
225             out->Vc[0]   = 0;
226
227             p_nbk_ener[coult][ljtreatment](nbl[nb], nbat,
228                                            ic,
229                                            shift_vec,
230                                            out->f,
231                                            fshift_p,
232                                            out->Vvdw,
233                                            out->Vc);
234         }}
235         else
236         {{
237             /* Calculate energy group contributions */
238             int i;
239
240             for (i = 0; i < out->nVS; i++)
241             {{
242                 out->VSvdw[i] = 0;
243             }}
244             for (i = 0; i < out->nVS; i++)
245             {{
246                 out->VSc[i] = 0;
247             }}
248
249             p_nbk_energrp[coult][ljtreatment](nbl[nb], nbat,
250                                               ic,
251                                               shift_vec,
252                                               out->f,
253                                               fshift_p,
254                                               out->VSvdw,
255                                               out->VSc);
256
257             reduce_group_energies(nbat->nenergrp, nbat->neg_2log,
258                                   out->VSvdw, out->VSc,
259                                   out->Vvdw, out->Vc);
260         }}
261     }}
262
263     if (force_flags & GMX_FORCE_ENERGY)
264     {{
265         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
266     }}
267 }}
268 #else
269 {{
270     gmx_incons("{5} called when such kernels "
271                " are not enabled.");
272 }}
273 #endif
274 #undef GMX_SIMD_J_UNROLL_SIZE