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 /* Analytical reaction-field kernels */
48 #define CALC_COUL_RF
49
50 /* Include the force+energy kernels */
51 #define CALC_ENERGIES
52 #include "nbnxn_kernel_ref_outer.h"
53 #undef CALC_ENERGIES
54
55 /* Include the force+energygroups kernels */
56 #define CALC_ENERGIES
57 #define ENERGY_GROUPS
58 #include "nbnxn_kernel_ref_outer.h"
59 #undef ENERGY_GROUPS
60 #undef CALC_ENERGIES
61
62 /* Include the force only kernels */
63 #include "nbnxn_kernel_ref_outer.h"
64
65 #undef CALC_COUL_RF
66
67
68 /* Tabulated exclusion interaction electrostatics kernels */
69 #define CALC_COUL_TAB
70
71 /* Include the force+energy kernels */
72 #define CALC_ENERGIES
73 #include "nbnxn_kernel_ref_outer.h"
74 #undef CALC_ENERGIES
75
76 /* Include the force+energygroups kernels */
77 #define CALC_ENERGIES
78 #define ENERGY_GROUPS
79 #include "nbnxn_kernel_ref_outer.h"
80 #undef ENERGY_GROUPS
81 #undef CALC_ENERGIES
82
83 /* Include the force only kernels */
84 #include "nbnxn_kernel_ref_outer.h"
85
86 /* Twin-range cut-off kernels */
87 #define VDW_CUTOFF_CHECK
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 #undef VDW_CUTOFF_CHECK
105
106 #undef CALC_COUL_TAB
107
108
109 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
110                                 const nbnxn_atomdata_t     *nbat,
111                                 const interaction_const_t  *ic,
112                                 rvec                       *shift_vec,
113                                 real                       *f,
114                                 real                       *fshift,
115                                 real                       *Vvdw,
116                                 real                       *Vc);
117
118 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t     *nbl,
119                                   const nbnxn_atomdata_t     *nbat,
120                                   const interaction_const_t  *ic,
121                                   rvec                       *shift_vec,
122                                   real                       *f,
123                                   real                       *fshift);
124
125 enum {
126     coultRF, coultTAB, coultTAB_TWIN, coultNR
127 };
128
129 p_nbk_func_ener p_nbk_c_ener[coultNR] =
130 {
131     nbnxn_kernel_ref_rf_ener,
132     nbnxn_kernel_ref_tab_ener,
133     nbnxn_kernel_ref_tab_twin_ener
134 };
135
136 p_nbk_func_ener p_nbk_c_energrp[coultNR] =
137 {
138     nbnxn_kernel_ref_rf_energrp,
139     nbnxn_kernel_ref_tab_energrp,
140     nbnxn_kernel_ref_tab_twin_energrp
141 };
142
143 p_nbk_func_noener p_nbk_c_noener[coultNR] =
144 {
145     nbnxn_kernel_ref_rf_noener,
146     nbnxn_kernel_ref_tab_noener,
147     nbnxn_kernel_ref_tab_twin_noener
148 };
149
150 void
151 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
152                  const nbnxn_atomdata_t     *nbat,
153                  const interaction_const_t  *ic,
154                  rvec                       *shift_vec,
155                  int                         force_flags,
156                  int                         clearF,
157                  real                       *fshift,
158                  real                       *Vc,
159                  real                       *Vvdw)
160 {
161     int                nnbl;
162     nbnxn_pairlist_t **nbl;
163     int                coult;
164     int                nb;
165
166     nnbl = nbl_list->nnbl;
167     nbl  = nbl_list->nbl;
168
169     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
170     {
171         coult = coultRF;
172     }
173     else
174     {
175         if (ic->rcoulomb == ic->rvdw)
176         {
177             coult = coultTAB;
178         }
179         else
180         {
181             coult = coultTAB_TWIN;
182         }
183     }
184
185 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
186     for (nb = 0; nb < nnbl; nb++)
187     {
188         nbnxn_atomdata_output_t *out;
189         real                    *fshift_p;
190
191         out = &nbat->out[nb];
192
193         if (clearF == enbvClearFYes)
194         {
195             clear_f(nbat, nb, out->f);
196         }
197
198         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
199         {
200             fshift_p = fshift;
201         }
202         else
203         {
204             fshift_p = out->fshift;
205
206             if (clearF == enbvClearFYes)
207             {
208                 clear_fshift(fshift_p);
209             }
210         }
211
212         if (!(force_flags & GMX_FORCE_ENERGY))
213         {
214             /* Don't calculate energies */
215             p_nbk_c_noener[coult](nbl[nb], nbat,
216                                   ic,
217                                   shift_vec,
218                                   out->f,
219                                   fshift_p);
220         }
221         else if (out->nV == 1)
222         {
223             /* No energy groups */
224             out->Vvdw[0] = 0;
225             out->Vc[0]   = 0;
226
227             p_nbk_c_ener[coult](nbl[nb], nbat,
228                                 ic,
229                                 shift_vec,
230                                 out->f,
231                                 fshift_p,
232                                 out->Vvdw,
233                                 out->Vc);
234         }
235         else
236         {
237             /* Calculate energy group contributions */
238             int i;
239
240             for (i = 0; i < out->nV; i++)
241             {
242                 out->Vvdw[i] = 0;
243             }
244             for (i = 0; i < out->nV; i++)
245             {
246                 out->Vc[i] = 0;
247             }
248
249             p_nbk_c_energrp[coult](nbl[nb], nbat,
250                                    ic,
251                                    shift_vec,
252                                    out->f,
253                                    fshift_p,
254                                    out->Vvdw,
255                                    out->Vc);
256         }
257     }
258
259     if (force_flags & GMX_FORCE_ENERGY)
260     {
261         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
262     }
263 }