Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_ref.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
9  * Copyright (c) 2001-2009, The GROMACS Development Team
10  *
11  * Gromacs is a library for molecular simulation and trajectory analysis,
12  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
13  * a full list of developers and information, check out http://www.gromacs.org
14  *
15  * This program is free software; you can redistribute it and/or modify it under
16  * the terms of the GNU Lesser General Public License as published by the Free
17  * Software Foundation; either version 2 of the License, or (at your option) any
18  * later version.
19  * As a special exception, you may use this file as part of a free software
20  * library without restriction.  Specifically, if other files instantiate
21  * templates or use macros or inline functions from this file, or you compile
22  * this file and link it with other files to produce an executable, this
23  * file does not by itself cause the resulting executable to be covered by
24  * the GNU Lesser General Public License.
25  *
26  * In plain-speak: do not worry about classes/macros/templates either - only
27  * changes to the library have to be LGPL, not an application linking with it.
28  *
29  * To help fund GROMACS development, we humbly ask that you cite
30  * the papers people have written on it - you can find them on the website!
31  */
32 #ifdef HAVE_CONFIG_H
33 #include <config.h>
34 #endif
35
36 #include <math.h>
37
38 #include "typedefs.h"
39 #include "vec.h"
40 #include "smalloc.h"
41 #include "force.h"
42 #include "gmx_omp_nthreads.h"
43 #include "nbnxn_kernel_ref.h"
44 #include "../nbnxn_consts.h"
45 #include "nbnxn_kernel_common.h"
46
47 /*! \brief Typedefs for declaring lookup tables of kernel functions.
48  */
49 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
50                                 const nbnxn_atomdata_t     *nbat,
51                                 const interaction_const_t  *ic,
52                                 rvec                       *shift_vec,
53                                 real                       *f,
54                                 real                       *fshift,
55                                 real                       *Vvdw,
56                                 real                       *Vc);
57
58 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t     *nbl,
59                                   const nbnxn_atomdata_t     *nbat,
60                                   const interaction_const_t  *ic,
61                                   rvec                       *shift_vec,
62                                   real                       *f,
63                                   real                       *fshift);
64
65 /* Analytical reaction-field kernels */
66 #define CALC_COUL_RF
67
68 /* Include the force+energy kernels */
69 #define CALC_ENERGIES
70 #include "nbnxn_kernel_ref_outer.h"
71 #undef CALC_ENERGIES
72
73 /* Include the force+energygroups kernels */
74 #define CALC_ENERGIES
75 #define ENERGY_GROUPS
76 #include "nbnxn_kernel_ref_outer.h"
77 #undef ENERGY_GROUPS
78 #undef CALC_ENERGIES
79
80 /* Include the force only kernels */
81 #include "nbnxn_kernel_ref_outer.h"
82
83 #undef CALC_COUL_RF
84
85
86 /* Tabulated exclusion interaction electrostatics kernels */
87 #define CALC_COUL_TAB
88
89 /* Include the force+energy kernels */
90 #define CALC_ENERGIES
91 #include "nbnxn_kernel_ref_outer.h"
92 #undef CALC_ENERGIES
93
94 /* Include the force+energygroups kernels */
95 #define CALC_ENERGIES
96 #define ENERGY_GROUPS
97 #include "nbnxn_kernel_ref_outer.h"
98 #undef ENERGY_GROUPS
99 #undef CALC_ENERGIES
100
101 /* Include the force only kernels */
102 #include "nbnxn_kernel_ref_outer.h"
103
104 /* Twin-range cut-off kernels */
105 #define VDW_CUTOFF_CHECK
106
107 /* Include the force+energy kernels */
108 #define CALC_ENERGIES
109 #include "nbnxn_kernel_ref_outer.h"
110 #undef CALC_ENERGIES
111
112 /* Include the force+energygroups kernels */
113 #define CALC_ENERGIES
114 #define ENERGY_GROUPS
115 #include "nbnxn_kernel_ref_outer.h"
116 #undef ENERGY_GROUPS
117 #undef CALC_ENERGIES
118
119 /* Include the force only kernels */
120 #include "nbnxn_kernel_ref_outer.h"
121
122 #undef VDW_CUTOFF_CHECK
123
124 #undef CALC_COUL_TAB
125
126
127 enum {
128     coultRF, coultTAB, coultTAB_TWIN, coultNR
129 };
130
131 p_nbk_func_ener p_nbk_c_ener[coultNR] =
132 {
133     nbnxn_kernel_ref_rf_ener,
134     nbnxn_kernel_ref_tab_ener,
135     nbnxn_kernel_ref_tab_twin_ener
136 };
137
138 p_nbk_func_ener p_nbk_c_energrp[coultNR] =
139 {
140     nbnxn_kernel_ref_rf_energrp,
141     nbnxn_kernel_ref_tab_energrp,
142     nbnxn_kernel_ref_tab_twin_energrp
143 };
144
145 p_nbk_func_noener p_nbk_c_noener[coultNR] =
146 {
147     nbnxn_kernel_ref_rf_noener,
148     nbnxn_kernel_ref_tab_noener,
149     nbnxn_kernel_ref_tab_twin_noener
150 };
151
152 void
153 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
154                  const nbnxn_atomdata_t     *nbat,
155                  const interaction_const_t  *ic,
156                  rvec                       *shift_vec,
157                  int                         force_flags,
158                  int                         clearF,
159                  real                       *fshift,
160                  real                       *Vc,
161                  real                       *Vvdw)
162 {
163     int                nnbl;
164     nbnxn_pairlist_t **nbl;
165     int                coult;
166     int                nb;
167
168     nnbl = nbl_list->nnbl;
169     nbl  = nbl_list->nbl;
170
171     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
172     {
173         coult = coultRF;
174     }
175     else
176     {
177         if (ic->rcoulomb == ic->rvdw)
178         {
179             coult = coultTAB;
180         }
181         else
182         {
183             coult = coultTAB_TWIN;
184         }
185     }
186
187 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
188     for (nb = 0; nb < nnbl; nb++)
189     {
190         nbnxn_atomdata_output_t *out;
191         real                    *fshift_p;
192
193         out = &nbat->out[nb];
194
195         if (clearF == enbvClearFYes)
196         {
197             clear_f(nbat, nb, out->f);
198         }
199
200         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
201         {
202             fshift_p = fshift;
203         }
204         else
205         {
206             fshift_p = out->fshift;
207
208             if (clearF == enbvClearFYes)
209             {
210                 clear_fshift(fshift_p);
211             }
212         }
213
214         if (!(force_flags & GMX_FORCE_ENERGY))
215         {
216             /* Don't calculate energies */
217             p_nbk_c_noener[coult](nbl[nb], nbat,
218                                   ic,
219                                   shift_vec,
220                                   out->f,
221                                   fshift_p);
222         }
223         else if (out->nV == 1)
224         {
225             /* No energy groups */
226             out->Vvdw[0] = 0;
227             out->Vc[0]   = 0;
228
229             p_nbk_c_ener[coult](nbl[nb], nbat,
230                                 ic,
231                                 shift_vec,
232                                 out->f,
233                                 fshift_p,
234                                 out->Vvdw,
235                                 out->Vc);
236         }
237         else
238         {
239             /* Calculate energy group contributions */
240             int i;
241
242             for (i = 0; i < out->nV; i++)
243             {
244                 out->Vvdw[i] = 0;
245             }
246             for (i = 0; i < out->nV; i++)
247             {
248                 out->Vc[i] = 0;
249             }
250
251             p_nbk_c_energrp[coult](nbl[nb], nbat,
252                                    ic,
253                                    shift_vec,
254                                    out->f,
255                                    fshift_p,
256                                    out->Vvdw,
257                                    out->Vc);
258         }
259     }
260
261     if (force_flags & GMX_FORCE_ENERGY)
262     {
263         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
264     }
265 }