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