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