Merge branch 'release-4-6'
[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, 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 "gmx_simd_macros.h"
45 #include "gmx_simd_vec.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
53 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
54  */
55 enum {{
56     coultRF, coultTAB, coultTAB_TWIN, coultEWALD, coultEWALD_TWIN, coultNR
57 }};
58
59 /* Declare and define the kernel function pointer lookup tables. */
60 static p_nbk_func_ener p_nbk_ener[coultNR][ljcrNR] =
61 {7}
62 static p_nbk_func_ener p_nbk_energrp[coultNR][ljcrNR] =
63 {8}
64 static p_nbk_func_noener p_nbk_noener[coultNR][ljcrNR] =
65 {9}
66
67 static void
68 reduce_group_energies(int ng, int ng_2log,
69                       const real *VSvdw, const real *VSc,
70                       real *Vvdw, real *Vc)
71 {{
72     const int unrollj      = GMX_SIMD_WIDTH_HERE/GMX_SIMD_J_UNROLL_SIZE;
73     const int unrollj_half = unrollj/2;
74     int       ng_p2, i, j, j0, j1, c, s;
75
76     ng_p2 = (1<<ng_2log);
77
78     /* The size of the x86 SIMD energy group buffer array is:
79      * ng*ng*ng_p2*unrollj_half*simd_width
80      */
81     for (i = 0; i < ng; i++)
82     {{
83         for (j = 0; j < ng; j++)
84         {{
85             Vvdw[i*ng+j] = 0;
86             Vc[i*ng+j]   = 0;
87         }}
88
89         for (j1 = 0; j1 < ng; j1++)
90         {{
91             for (j0 = 0; j0 < ng; j0++)
92             {{
93                 c = ((i*ng + j1)*ng_p2 + j0)*unrollj_half*unrollj;
94                 for (s = 0; s < unrollj_half; s++)
95                 {{
96                     Vvdw[i*ng+j0] += VSvdw[c+0];
97                     Vvdw[i*ng+j1] += VSvdw[c+1];
98                     Vc  [i*ng+j0] += VSc  [c+0];
99                     Vc  [i*ng+j1] += VSc  [c+1];
100                     c             += unrollj + 2;
101                 }}
102             }}
103         }}
104     }}
105 }}
106
107 #else /* {0} */
108
109 #include "gmx_fatal.h"
110
111 #endif /* {0} */
112
113 void
114 {5}(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
115 {6}const nbnxn_atomdata_t    gmx_unused *nbat,
116 {6}const interaction_const_t gmx_unused *ic,
117 {6}int                       gmx_unused  ewald_excl,
118 {6}rvec                      gmx_unused *shift_vec,
119 {6}int                       gmx_unused  force_flags,
120 {6}int                       gmx_unused  clearF,
121 {6}real                      gmx_unused *fshift,
122 {6}real                      gmx_unused *Vc,
123 {6}real                      gmx_unused *Vvdw)
124 #ifdef {0}
125 {{
126     int                nnbl;
127     nbnxn_pairlist_t **nbl;
128     int                coult;
129     int                nb;
130
131     nnbl = nbl_list->nnbl;
132     nbl  = nbl_list->nbl;
133
134     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
135     {{
136         coult = coultRF;
137     }}
138     else
139     {{
140         if (ewald_excl == ewaldexclTable)
141         {{
142             if (ic->rcoulomb == ic->rvdw)
143             {{
144                 coult = coultTAB;
145             }}
146             else
147             {{
148                 coult = coultTAB_TWIN;
149             }}
150         }}
151         else
152         {{
153             if (ic->rcoulomb == ic->rvdw)
154             {{
155                 coult = coultEWALD;
156             }}
157             else
158             {{
159                 coult = coultEWALD_TWIN;
160             }}
161         }}
162     }}
163
164 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
165     for (nb = 0; nb < nnbl; nb++)
166     {{
167         nbnxn_atomdata_output_t *out;
168         real                    *fshift_p;
169
170         out = &nbat->out[nb];
171
172         if (clearF == enbvClearFYes)
173         {{
174             clear_f(nbat, nb, out->f);
175         }}
176
177         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
178         {{
179             fshift_p = fshift;
180         }}
181         else
182         {{
183             fshift_p = out->fshift;
184
185             if (clearF == enbvClearFYes)
186             {{
187                 clear_fshift(fshift_p);
188             }}
189         }}
190
191         if (!(force_flags & GMX_FORCE_ENERGY))
192         {{
193             /* Don't calculate energies */
194             p_nbk_noener[coult][nbat->comb_rule](nbl[nb], nbat,
195                                                  ic,
196                                                  shift_vec,
197                                                  out->f,
198                                                  fshift_p);
199         }}
200         else if (out->nV == 1)
201         {{
202             /* No energy groups */
203             out->Vvdw[0] = 0;
204             out->Vc[0]   = 0;
205
206             p_nbk_ener[coult][nbat->comb_rule](nbl[nb], nbat,
207                                                ic,
208                                                shift_vec,
209                                                out->f,
210                                                fshift_p,
211                                                out->Vvdw,
212                                                out->Vc);
213         }}
214         else
215         {{
216             /* Calculate energy group contributions */
217             int i;
218
219             for (i = 0; i < out->nVS; i++)
220             {{
221                 out->VSvdw[i] = 0;
222             }}
223             for (i = 0; i < out->nVS; i++)
224             {{
225                 out->VSc[i] = 0;
226             }}
227
228             p_nbk_energrp[coult][nbat->comb_rule](nbl[nb], nbat,
229                                                   ic,
230                                                   shift_vec,
231                                                   out->f,
232                                                   fshift_p,
233                                                   out->VSvdw,
234                                                   out->VSc);
235
236             reduce_group_energies(nbat->nenergrp, nbat->neg_2log,
237                                   out->VSvdw, out->VSc,
238                                   out->Vvdw, out->Vc);
239         }}
240     }}
241
242     if (force_flags & GMX_FORCE_ENERGY)
243     {{
244         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
245     }}
246 }}
247 #else
248 {{
249     gmx_incons("{5} called when such kernels "
250                " are not enabled.");
251 }}
252 #endif
253 #undef GMX_SIMD_J_UNROLL_SIZE