Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / nbnxm / kernel_common.h
blob7a2b97a7aa366528386fdec7e215833eea7fcc41
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017,2018,2019, 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 nbnxm pair interaction kernel function types and kind counts, also declares utility functions used in nbnxm_kernel.cpp.
41 * \author Berk Hess <hess@kth.se>
44 #ifndef GMX_NBXNM_KERNEL_COMMON_H
45 #define GMX_NBXNM_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/nbnxm/atomdata.h"
50 #include "gromacs/nbnxm/pairlist.h"
51 #include "gromacs/utility/real.h"
53 struct interaction_const_t;
55 // TODO: Consider using one nbk_func type now ener and noener are identical
57 /*! \brief Pair-interaction kernel type that also calculates energies.
59 typedef void (nbk_func_ener)(const NbnxnPairlistCpu *nbl,
60 const nbnxn_atomdata_t *nbat,
61 const interaction_const_t *ic,
62 const rvec *shift_vec,
63 nbnxn_atomdata_output_t *out);
65 /*! \brief Pointer to \p nbk_func_ener.
67 typedef nbk_func_ener *p_nbk_func_ener;
69 /*! \brief Pair-interaction kernel type that does not calculates energies.
71 typedef void (nbk_func_noener)(const NbnxnPairlistCpu *nbl,
72 const nbnxn_atomdata_t *nbat,
73 const interaction_const_t *ic,
74 const rvec *shift_vec,
75 nbnxn_atomdata_output_t *out);
77 /*! \brief Pointer to \p nbk_func_noener.
79 typedef nbk_func_noener *p_nbk_func_noener;
81 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
83 enum {
84 coulktRF, coulktTAB, coulktTAB_TWIN, coulktEWALD, coulktEWALD_TWIN, coulktNR
87 /*! \brief Kinds of Van der Waals treatments in SIMD Verlet kernels
89 * The \p LJCUT_COMB refers to the LJ combination rule for the short range.
90 * The \p EWALDCOMB refers to the combination rule for the grid part.
91 * \p vdwktNR is the number of VdW treatments for the SIMD kernels.
92 * \p vdwktNR_ref is the number of VdW treatments for the C reference kernels.
93 * These two numbers differ, because currently only the reference kernels
94 * support LB combination rules for the LJ-Ewald grid part.
96 enum {
97 vdwktLJCUT_COMBGEOM, vdwktLJCUT_COMBLB, vdwktLJCUT_COMBNONE, vdwktLJFORCESWITCH, vdwktLJPOTSWITCH, vdwktLJEWALDCOMBGEOM, vdwktLJEWALDCOMBLB, vdwktNR = vdwktLJEWALDCOMBLB, vdwktNR_ref
100 /*! \brief Clears the force buffer.
102 * Either the whole buffer is cleared or only the parts used
103 * by thread/task \p outputIndex when nbat->bUseBufferFlags is set.
105 * \param[in,out] nbat The Nbnxm atom data
106 * \param[in] outputIndex The index of the output object to clear
108 void
109 clearForceBuffer(nbnxn_atomdata_t *nbat,
110 int outputIndex);
112 /*! \brief Clears the shift forces.
114 void
115 clear_fshift(real *fshift);
117 /*! \brief Reduces the collected energy terms over the pair-lists/threads.
119 void
120 reduce_energies_over_lists(const nbnxn_atomdata_t *nbat,
121 int nlist,
122 real *Vvdw,
123 real *Vc);
125 #endif