efcafd3f1d857ecb2325586405780198d63bde74
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_ref.c
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 "gmxpre.h"
36
37 #include "config.h"
38
39 #include <math.h>
40 #include <assert.h>
41
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/legacyheaders/force.h"
44 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
45
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/ishift.h"
48 #include "gromacs/utility/smalloc.h"
49
50 #include "nbnxn_kernel_ref.h"
51 #include "../nbnxn_consts.h"
52 #include "nbnxn_kernel_common.h"
53 #include "gromacs/mdlib/nb_verlet.h"
54
55 /*! \brief Typedefs for declaring lookup tables of kernel functions.
56  */
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 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t     *nbl,
66                                 const nbnxn_atomdata_t     *nbat,
67                                 const interaction_const_t  *ic,
68                                 rvec                       *shift_vec,
69                                 real                       *f,
70                                 real                       *fshift,
71                                 real                       *Vvdw,
72                                 real                       *Vc);
73
74 /* Analytical reaction-field kernels */
75 #define CALC_COUL_RF
76 #define LJ_CUT
77 #include "nbnxn_kernel_ref_includes.h"
78 #undef LJ_CUT
79 #define LJ_FORCE_SWITCH
80 #include "nbnxn_kernel_ref_includes.h"
81 #undef LJ_FORCE_SWITCH
82 #define LJ_POT_SWITCH
83 #include "nbnxn_kernel_ref_includes.h"
84 #undef LJ_POT_SWITCH
85 #define LJ_EWALD
86 #define LJ_CUT
87 #define LJ_EWALD_COMB_GEOM
88 #include "nbnxn_kernel_ref_includes.h"
89 #undef LJ_EWALD_COMB_GEOM
90 #define LJ_EWALD_COMB_LB
91 #include "nbnxn_kernel_ref_includes.h"
92 #undef LJ_EWALD_COMB_LB
93 #undef LJ_CUT
94 #undef LJ_EWALD
95 #undef CALC_COUL_RF
96
97
98 /* Tabulated exclusion interaction electrostatics kernels */
99 #define CALC_COUL_TAB
100 #define LJ_CUT
101 #include "nbnxn_kernel_ref_includes.h"
102 #undef LJ_CUT
103 #define LJ_FORCE_SWITCH
104 #include "nbnxn_kernel_ref_includes.h"
105 #undef LJ_FORCE_SWITCH
106 #define LJ_POT_SWITCH
107 #include "nbnxn_kernel_ref_includes.h"
108 #undef LJ_POT_SWITCH
109 #define LJ_EWALD
110 #define LJ_CUT
111 #define LJ_EWALD_COMB_GEOM
112 #include "nbnxn_kernel_ref_includes.h"
113 #undef LJ_EWALD_COMB_GEOM
114 #define LJ_EWALD_COMB_LB
115 #include "nbnxn_kernel_ref_includes.h"
116 #undef LJ_EWALD_COMB_LB
117 #undef LJ_CUT
118 #undef LJ_EWALD
119 /* Twin-range cut-off kernels */
120 #define VDW_CUTOFF_CHECK
121 #define LJ_CUT
122 #include "nbnxn_kernel_ref_includes.h"
123 #undef LJ_CUT
124 #define LJ_FORCE_SWITCH
125 #include "nbnxn_kernel_ref_includes.h"
126 #undef LJ_FORCE_SWITCH
127 #define LJ_POT_SWITCH
128 #include "nbnxn_kernel_ref_includes.h"
129 #undef LJ_POT_SWITCH
130 #define LJ_EWALD
131 #define LJ_CUT
132 #define LJ_EWALD_COMB_GEOM
133 #include "nbnxn_kernel_ref_includes.h"
134 #undef LJ_EWALD_COMB_GEOM
135 #define LJ_EWALD_COMB_LB
136 #include "nbnxn_kernel_ref_includes.h"
137 #undef LJ_EWALD_COMB_LB
138 #undef LJ_CUT
139 #undef LJ_EWALD
140 #undef VDW_CUTOFF_CHECK
141 #undef CALC_COUL_TAB
142
143
144 enum {
145     coultRF, coultTAB, coultTAB_TWIN, coultNR
146 };
147
148 enum {
149     vdwtCUT, vdwtFSWITCH, vdwtPSWITCH, vdwtEWALDGEOM, vdwtEWALDLB, vdwtNR
150 };
151
152 p_nbk_func_noener p_nbk_c_noener[coultNR][vdwtNR] =
153 {
154     { nbnxn_kernel_ElecRF_VdwLJ_F_ref,           nbnxn_kernel_ElecRF_VdwLJFsw_F_ref,           nbnxn_kernel_ElecRF_VdwLJPsw_F_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_ref           },
155     { nbnxn_kernel_ElecQSTab_VdwLJ_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJFsw_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJPsw_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_F_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombLB_F_ref        },
156     { nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJFsw_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJPsw_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_F_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombLB_F_ref }
157 };
158
159 p_nbk_func_ener p_nbk_c_ener[coultNR][vdwtNR] =
160 {
161     { nbnxn_kernel_ElecRF_VdwLJ_VF_ref,           nbnxn_kernel_ElecRF_VdwLJFsw_VF_ref,           nbnxn_kernel_ElecRF_VdwLJPsw_VF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_ref            },
162     { nbnxn_kernel_ElecQSTab_VdwLJ_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJFsw_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJPsw_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombLB_VF_ref         },
163     { nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJFsw_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJPsw_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombLB_VF_ref  }
164 };
165
166 p_nbk_func_ener p_nbk_c_energrp[coultNR][vdwtNR] =
167 {
168     { nbnxn_kernel_ElecRF_VdwLJ_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJFsw_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJPsw_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VgrpF_ref,           nbnxn_kernel_ElecRF_VdwLJEwCombLB_VgrpF_ref           },
169     { nbnxn_kernel_ElecQSTab_VdwLJ_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJFsw_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJPsw_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VgrpF_ref,        nbnxn_kernel_ElecQSTab_VdwLJEwCombLB_VgrpF_ref        },
170     { nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJFsw_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJPsw_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VgrpF_ref, nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombLB_VgrpF_ref }
171 };
172
173 void
174 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
175                  const nbnxn_atomdata_t     *nbat,
176                  const interaction_const_t  *ic,
177                  rvec                       *shift_vec,
178                  int                         force_flags,
179                  int                         clearF,
180                  real                       *fshift,
181                  real                       *Vc,
182                  real                       *Vvdw)
183 {
184     int                nnbl;
185     nbnxn_pairlist_t **nbl;
186     int                coult;
187     int                vdwt;
188     int                nb;
189
190     nnbl = nbl_list->nnbl;
191     nbl  = nbl_list->nbl;
192
193     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
194     {
195         coult = coultRF;
196     }
197     else
198     {
199         if (ic->rcoulomb == ic->rvdw)
200         {
201             coult = coultTAB;
202         }
203         else
204         {
205             coult = coultTAB_TWIN;
206         }
207     }
208
209     if (ic->vdwtype == evdwCUT)
210     {
211         switch (ic->vdw_modifier)
212         {
213             case eintmodPOTSHIFT:
214             case eintmodNONE:
215                 vdwt = vdwtCUT;
216                 break;
217             case eintmodFORCESWITCH:
218                 vdwt = vdwtFSWITCH;
219                 break;
220             case eintmodPOTSWITCH:
221                 vdwt = vdwtPSWITCH;
222                 break;
223             default:
224                 gmx_incons("Unsupported VdW modifier");
225                 break;
226         }
227     }
228     else if (ic->vdwtype == evdwPME)
229     {
230         if (ic->ljpme_comb_rule == ljcrGEOM)
231         {
232             assert(nbat->comb_rule == ljcrGEOM);
233             vdwt = vdwtEWALDGEOM;
234         }
235         else
236         {
237             assert(nbat->comb_rule == ljcrLB);
238             vdwt = vdwtEWALDLB;
239         }
240     }
241     else
242     {
243         gmx_incons("Unsupported vdwtype in nbnxn reference kernel");
244     }
245
246 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
247     for (nb = 0; nb < nnbl; nb++)
248     {
249         nbnxn_atomdata_output_t *out;
250         real                    *fshift_p;
251
252         out = &nbat->out[nb];
253
254         if (clearF == enbvClearFYes)
255         {
256             clear_f(nbat, nb, out->f);
257         }
258
259         if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
260         {
261             fshift_p = fshift;
262         }
263         else
264         {
265             fshift_p = out->fshift;
266
267             if (clearF == enbvClearFYes)
268             {
269                 clear_fshift(fshift_p);
270             }
271         }
272
273         if (!(force_flags & GMX_FORCE_ENERGY))
274         {
275             /* Don't calculate energies */
276             p_nbk_c_noener[coult][vdwt](nbl[nb], nbat,
277                                         ic,
278                                         shift_vec,
279                                         out->f,
280                                         fshift_p);
281         }
282         else if (out->nV == 1)
283         {
284             /* No energy groups */
285             out->Vvdw[0] = 0;
286             out->Vc[0]   = 0;
287
288             p_nbk_c_ener[coult][vdwt](nbl[nb], nbat,
289                                       ic,
290                                       shift_vec,
291                                       out->f,
292                                       fshift_p,
293                                       out->Vvdw,
294                                       out->Vc);
295         }
296         else
297         {
298             /* Calculate energy group contributions */
299             int i;
300
301             for (i = 0; i < out->nV; i++)
302             {
303                 out->Vvdw[i] = 0;
304             }
305             for (i = 0; i < out->nV; i++)
306             {
307                 out->Vc[i] = 0;
308             }
309
310             p_nbk_c_energrp[coult][vdwt](nbl[nb], nbat,
311                                          ic,
312                                          shift_vec,
313                                          out->f,
314                                          fshift_p,
315                                          out->Vvdw,
316                                          out->Vc);
317         }
318     }
319
320     if (force_flags & GMX_FORCE_ENERGY)
321     {
322         reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
323     }
324 }