Unify nbnxn kernel dispatchers
[gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_common.h
blob4edc4c2495684c1c87166ddafba07127f522b53d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017, 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.
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.
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.
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.
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.
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.
36 /*! \internal \file
38 * \brief
39 * Declares the nbnxn pair interaction kernel function types and kind counts, also declares utility functions used in nbnxn_kernel.cpp.
41 * \author Berk Hess <hess@kth.se>
44 #ifndef _nbnxn_kernel_common_h
45 #define _nbnxn_kernel_common_h
47 #include "gromacs/math/vectypes.h"
48 /* nbnxn_atomdata_t and nbnxn_pairlist_t could be forward declared, but that requires modifications in all SIMD kernel files */
49 #include "gromacs/mdlib/nbnxn_atomdata.h"
50 #include "gromacs/mdlib/nbnxn_pairlist.h"
51 #include "gromacs/utility/real.h"
53 struct interaction_const_t;
55 /*! \brief Pair-interaction kernel type that also calculates energies.
57 typedef void (nbk_func_ener)(const nbnxn_pairlist_t *nbl,
58 const nbnxn_atomdata_t *nbat,
59 const interaction_const_t *ic,
60 rvec *shift_vec,
61 real *f,
62 real *fshift,
63 real *Vvdw,
64 real *Vc);
66 /*! \brief Pointer to \p nbk_func_ener.
68 typedef nbk_func_ener *p_nbk_func_ener;
70 /*! \brief Pair-interaction kernel type that does not calculates energies.
72 typedef void (nbk_func_noener)(const nbnxn_pairlist_t *nbl,
73 const nbnxn_atomdata_t *nbat,
74 const interaction_const_t *ic,
75 rvec *shift_vec,
76 real *f,
77 real *fshift);
79 /*! \brief Pointer to \p nbk_func_noener.
81 typedef nbk_func_noener *p_nbk_func_noener;
83 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
85 enum {
86 coulktRF, coulktTAB, coulktTAB_TWIN, coulktEWALD, coulktEWALD_TWIN, coulktNR
89 /*! \brief Kinds of Van der Waals treatments in SIMD Verlet kernels
91 * The \p LJCUT_COMB refers to the LJ combination rule for the short range.
92 * The \p EWALDCOMB refers to the combination rule for the grid part.
93 * \p vdwktNR is the number of VdW treatments for the SIMD kernels.
94 * \p vdwktNR_ref is the number of VdW treatments for the C reference kernels.
95 * These two numbers differ, because currently only the reference kernels
96 * support LB combination rules for the LJ-Ewald grid part.
98 enum {
99 vdwktLJCUT_COMBGEOM, vdwktLJCUT_COMBLB, vdwktLJCUT_COMBNONE, vdwktLJFORCESWITCH, vdwktLJPOTSWITCH, vdwktLJEWALDCOMBGEOM, vdwktLJEWALDCOMBLB, vdwktNR = vdwktLJEWALDCOMBLB, vdwktNR_ref
102 /*! \brief Clears the force buffer.
104 * Either the whole buffer is cleared or only the parts used
105 * by the current thread when nbat->bUseBufferFlags is set.
106 * In the latter case output_index is the task/thread list/buffer index.
108 void
109 clear_f(const nbnxn_atomdata_t *nbat, int output_index, real *f);
111 /*! \brief Clears the shift forces.
113 void
114 clear_fshift(real *fshift);
116 /*! \brief Reduces the collected energy terms over the pair-lists/threads.
118 void
119 reduce_energies_over_lists(const nbnxn_atomdata_t *nbat,
120 int nlist,
121 real *Vvdw,
122 real *Vc);
124 #endif