Implemented LJ-PME nbnxn kernels
[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     coulktRF, coulktTAB, coulktTAB_TWIN, coulktEWALD, coulktEWALD_TWIN, coulktNR
58 }};
59
60 /*! \brief Kinds of Van der Waals treatments in SIMD Verlet kernels
61  */
62 enum {{
63     vdwktLJCUT_COMBGEOM, vdwktLJCUT_COMBLB, vdwktLJCUT_COMBNONE, vdwktLJFORCESWITCH, vdwktLJPOTSWITCH, vdwktLJEWALDCOMBGEOM, vdwktNR
64 }};
65
66 /* Declare and define the kernel function pointer lookup tables.
67  * The minor index of the array goes over both the LJ combination rules,
68  * which is only supported by plain cut-off, and the LJ switch/PME functions.
69  */
70 static p_nbk_func_noener p_nbk_noener[coulktNR][vdwktNR] =
71 {7}
72 static p_nbk_func_ener p_nbk_ener[coulktNR][vdwktNR] =
73 {8}
74 static p_nbk_func_ener p_nbk_energrp[coulktNR][vdwktNR] =
75 {9}
76
77 static void
78 reduce_group_energies(int ng, int ng_2log,
79                       const real *VSvdw, const real *VSc,
80                       real *Vvdw, real *Vc)
81 {{
82     const int unrollj      = GMX_SIMD_REAL_WIDTH/GMX_SIMD_J_UNROLL_SIZE;
83     const int unrollj_half = unrollj/2;
84     int       ng_p2, i, j, j0, j1, c, s;
85
86     ng_p2 = (1<<ng_2log);
87
88     /* The size of the x86 SIMD energy group buffer array is:
89      * ng*ng*ng_p2*unrollj_half*simd_width
90      */
91     for (i = 0; i < ng; i++)
92     {{
93         for (j = 0; j < ng; j++)
94         {{
95             Vvdw[i*ng+j] = 0;
96             Vc[i*ng+j]   = 0;
97         }}
98
99         for (j1 = 0; j1 < ng; j1++)
100         {{
101             for (j0 = 0; j0 < ng; j0++)
102             {{
103                 c = ((i*ng + j1)*ng_p2 + j0)*unrollj_half*unrollj;
104                 for (s = 0; s < unrollj_half; s++)
105                 {{
106                     Vvdw[i*ng+j0] += VSvdw[c+0];
107                     Vvdw[i*ng+j1] += VSvdw[c+1];
108                     Vc  [i*ng+j0] += VSc  [c+0];
109                     Vc  [i*ng+j1] += VSc  [c+1];
110                     c             += unrollj + 2;
111                 }}
112             }}
113         }}
114     }}
115 }}
116
117 #else /* {0} */
118
119 #include "gmx_fatal.h"
120
121 #endif /* {0} */
122
123 void
124 {5}(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
125 {6}const nbnxn_atomdata_t    gmx_unused *nbat,
126 {6}const interaction_const_t gmx_unused *ic,
127 {6}int                       gmx_unused  ewald_excl,
128 {6}rvec                      gmx_unused *shift_vec,
129 {6}int                       gmx_unused  force_flags,
130 {6}int                       gmx_unused  clearF,
131 {6}real                      gmx_unused *fshift,
132 {6}real                      gmx_unused *Vc,
133 {6}real                      gmx_unused *Vvdw)
134 #ifdef {0}
135 {{
136     int                nnbl;
137     nbnxn_pairlist_t **nbl;
138     int                coulkt, vdwkt = 0;
139     int                nb;
140
141     nnbl = nbl_list->nnbl;
142     nbl  = nbl_list->nbl;
143
144     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
145     {{
146         coulkt = coulktRF;
147     }}
148     else
149     {{
150         if (ewald_excl == ewaldexclTable)
151         {{
152             if (ic->rcoulomb == ic->rvdw)
153             {{
154                 coulkt = coulktTAB;
155             }}
156             else
157             {{
158                 coulkt = coulktTAB_TWIN;
159             }}
160         }}
161         else
162         {{
163             if (ic->rcoulomb == ic->rvdw)
164             {{
165                 coulkt = coulktEWALD;
166             }}
167             else
168             {{
169                 coulkt = coulktEWALD_TWIN;
170             }}
171         }}
172     }}
173
174     if (ic->vdwtype == evdwCUT)
175     {{
176         switch (ic->vdw_modifier)
177         {{
178             case eintmodNONE:
179             case eintmodPOTSHIFT:
180                 switch (nbat->comb_rule)
181                 {{
182                     case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
183                     case ljcrLB:   vdwkt = vdwktLJCUT_COMBLB;   break;
184                     case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
185                     default:       gmx_incons("Unknown combination rule");
186                 }}
187                 break;
188             case eintmodFORCESWITCH:
189                 vdwkt = vdwktLJFORCESWITCH;
190                 break;
191             case eintmodPOTSWITCH:
192                 vdwkt = vdwktLJPOTSWITCH;
193                 break;
194             default:
195                 gmx_incons("Unsupported VdW interaction modifier");
196         }}
197     }}
198     else if (ic->vdwtype == evdwPME)
199     {{
200         if (ic->ljpme_comb_rule == eljpmeLB)
201         {{
202             gmx_incons("The nbnxn SIMD kernels don't suport LJ-PME with LB");
203         }}
204         vdwkt = vdwktLJEWALDCOMBGEOM;
205     }}
206     else
207     {{
208         gmx_incons("Unsupported VdW interaction type");
209     }}
210
211 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
212     for (nb = 0; nb < nnbl; nb++)
213     {{
214         nbnxn_atomdata_output_t *out;
215         real                    *fshift_p;
216
217         out = &nbat->out[nb];
218
219         if (clearF == enbvClearFYes)
220         {{
221             clear_f(nbat, nb, out->f);
222         }}
223
224         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
225         {{
226             fshift_p = fshift;
227         }}
228         else
229         {{
230             fshift_p = out->fshift;
231
232             if (clearF == enbvClearFYes)
233             {{
234                 clear_fshift(fshift_p);
235             }}
236         }}
237
238         if (!(force_flags & GMX_FORCE_ENERGY))
239         {{
240             /* Don't calculate energies */
241             p_nbk_noener[coulkt][vdwkt](nbl[nb], nbat,
242                                         ic,
243                                         shift_vec,
244                                         out->f,
245                                         fshift_p);
246         }}
247         else if (out->nV == 1)
248         {{
249             /* No energy groups */
250             out->Vvdw[0] = 0;
251             out->Vc[0]   = 0;
252
253             p_nbk_ener[coulkt][vdwkt](nbl[nb], nbat,
254                                       ic,
255                                       shift_vec,
256                                       out->f,
257                                       fshift_p,
258                                       out->Vvdw,
259                                       out->Vc);
260         }}
261         else
262         {{
263             /* Calculate energy group contributions */
264             int i;
265
266             for (i = 0; i < out->nVS; i++)
267             {{
268                 out->VSvdw[i] = 0;
269             }}
270             for (i = 0; i < out->nVS; i++)
271             {{
272                 out->VSc[i] = 0;
273             }}
274
275             p_nbk_energrp[coulkt][vdwkt](nbl[nb], nbat,
276                                          ic,
277                                          shift_vec,
278                                          out->f,
279                                          fshift_p,
280                                          out->VSvdw,
281                                          out->VSc);
282
283             reduce_group_energies(nbat->nenergrp, nbat->neg_2log,
284                                   out->VSvdw, out->VSc,
285                                   out->Vvdw, out->Vc);
286         }}
287     }}
288
289     if (force_flags & GMX_FORCE_ENERGY)
290     {{
291         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
292     }}
293 }}
294 #else
295 {{
296     gmx_incons("{5} called when such kernels "
297                " are not enabled.");
298 }}
299 #endif
300 #undef GMX_SIMD_J_UNROLL_SIZE