added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / src / mdlib / nbnxn_kernels / nbnxn_kernel_ref.c
blobc1f9e2e40a8080b4e324188b939c47a60ae77a5d
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
6 * G R O M A C S
8 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
9 * Copyright (c) 2001-2009, The GROMACS Development Team
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
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.
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.
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!
32 #ifdef HAVE_CONFIG_H
33 #include <config.h>
34 #endif
36 #include <math.h>
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"
47 /* Analytical reaction-field kernels */
48 #define CALC_COUL_RF
50 /* Include the force+energy kernels */
51 #define CALC_ENERGIES
52 #include "nbnxn_kernel_ref_outer.h"
53 #undef CALC_ENERGIES
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
62 /* Include the force only kernels */
63 #include "nbnxn_kernel_ref_outer.h"
65 #undef CALC_COUL_RF
68 /* Tabulated exclusion interaction electrostatics kernels */
69 #define CALC_COUL_TAB
71 /* Include the force+energy kernels */
72 #define CALC_ENERGIES
73 #include "nbnxn_kernel_ref_outer.h"
74 #undef CALC_ENERGIES
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
83 /* Include the force only kernels */
84 #include "nbnxn_kernel_ref_outer.h"
86 /* Twin-range cut-off kernels */
87 #define VDW_CUTOFF_CHECK
89 /* Include the force+energy kernels */
90 #define CALC_ENERGIES
91 #include "nbnxn_kernel_ref_outer.h"
92 #undef CALC_ENERGIES
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
101 /* Include the force only kernels */
102 #include "nbnxn_kernel_ref_outer.h"
104 #undef VDW_CUTOFF_CHECK
106 #undef CALC_COUL_TAB
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);
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);
125 enum { coultRF, coultTAB, coultTAB_TWIN, coultNR };
127 p_nbk_func_ener p_nbk_c_ener[coultNR] =
128 { nbnxn_kernel_ref_rf_ener,
129 nbnxn_kernel_ref_tab_ener,
130 nbnxn_kernel_ref_tab_twin_ener };
132 p_nbk_func_ener p_nbk_c_energrp[coultNR] =
133 { nbnxn_kernel_ref_rf_energrp,
134 nbnxn_kernel_ref_tab_energrp,
135 nbnxn_kernel_ref_tab_twin_energrp};
137 p_nbk_func_noener p_nbk_c_noener[coultNR] =
138 { nbnxn_kernel_ref_rf_noener,
139 nbnxn_kernel_ref_tab_noener,
140 nbnxn_kernel_ref_tab_twin_noener };
142 void
143 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
144 const nbnxn_atomdata_t *nbat,
145 const interaction_const_t *ic,
146 rvec *shift_vec,
147 int force_flags,
148 int clearF,
149 real *fshift,
150 real *Vc,
151 real *Vvdw)
153 int nnbl;
154 nbnxn_pairlist_t **nbl;
155 int coult;
156 int nb;
158 nnbl = nbl_list->nnbl;
159 nbl = nbl_list->nbl;
161 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
163 coult = coultRF;
165 else
167 if (ic->rcoulomb == ic->rvdw)
169 coult = coultTAB;
171 else
173 coult = coultTAB_TWIN;
177 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
178 for(nb=0; nb<nnbl; nb++)
180 nbnxn_atomdata_output_t *out;
181 real *fshift_p;
183 out = &nbat->out[nb];
185 if (clearF == enbvClearFYes)
187 clear_f(nbat,out->f);
190 if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
192 fshift_p = fshift;
194 else
196 fshift_p = out->fshift;
198 if (clearF == enbvClearFYes)
200 clear_fshift(fshift_p);
204 if (!(force_flags & GMX_FORCE_ENERGY))
206 /* Don't calculate energies */
207 p_nbk_c_noener[coult](nbl[nb],nbat,
209 shift_vec,
210 out->f,
211 fshift_p);
213 else if (out->nV == 1)
215 /* No energy groups */
216 out->Vvdw[0] = 0;
217 out->Vc[0] = 0;
219 p_nbk_c_ener[coult](nbl[nb],nbat,
221 shift_vec,
222 out->f,
223 fshift_p,
224 out->Vvdw,
225 out->Vc);
227 else
229 /* Calculate energy group contributions */
230 int i;
232 for(i=0; i<out->nV; i++)
234 out->Vvdw[i] = 0;
236 for(i=0; i<out->nV; i++)
238 out->Vc[i] = 0;
241 p_nbk_c_energrp[coult](nbl[nb],nbat,
243 shift_vec,
244 out->f,
245 fshift_p,
246 out->Vvdw,
247 out->Vc);
251 if (force_flags & GMX_FORCE_ENERGY)
253 /* Reduce the energies */
254 for(nb=0; nb<nnbl; nb++)
256 int i;
258 for(i=0; i<nbat->out[nb].nV; i++)
260 Vvdw[i] += nbat->out[nb].Vvdw[i];
261 Vc[i] += nbat->out[nb].Vc[i];