f107ea591836a77d5789886622fa0af9280251a5
[alexxy/gromacs.git] / src / gromacs / nbnxm / kernel_common.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *
39  * \brief
40  * Declares the nbnxm pair interaction kernel function types and kind counts, also declares utility functions used in nbnxm_kernel.cpp.
41  *
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_nbnxm
44  */
45
46 #ifndef GMX_NBXNM_KERNEL_COMMON_H
47 #define GMX_NBXNM_KERNEL_COMMON_H
48
49 #include "gromacs/math/vectypes.h"
50 /* nbnxn_atomdata_t and nbnxn_pairlist_t could be forward declared, but that requires modifications in all SIMD kernel files */
51 #include "gromacs/nbnxm/atomdata.h"
52 #include "gromacs/utility/real.h"
53
54 #include "pairlist.h"
55
56 struct interaction_const_t;
57 enum class CoulombInteractionType : int;
58
59 namespace Nbnxm
60 {
61 enum class EwaldExclusionType : int;
62 }
63
64 // TODO: Consider using one nbk_func type now ener and noener are identical
65
66 /*! \brief Pair-interaction kernel type that also calculates energies.
67  */
68 typedef void(nbk_func_ener)(const NbnxnPairlistCpu*    nbl,
69                             const nbnxn_atomdata_t*    nbat,
70                             const interaction_const_t* ic,
71                             const rvec*                shift_vec,
72                             nbnxn_atomdata_output_t*   out);
73
74 /*! \brief Pointer to \p nbk_func_ener.
75  */
76 typedef nbk_func_ener* p_nbk_func_ener;
77
78 /*! \brief Pair-interaction kernel type that does not calculates energies.
79  */
80 typedef void(nbk_func_noener)(const NbnxnPairlistCpu*    nbl,
81                               const nbnxn_atomdata_t*    nbat,
82                               const interaction_const_t* ic,
83                               const rvec*                shift_vec,
84                               nbnxn_atomdata_output_t*   out);
85
86 /*! \brief Pointer to \p nbk_func_noener.
87  */
88 typedef nbk_func_noener* p_nbk_func_noener;
89
90 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
91  */
92 enum class CoulombKernelType : int
93 {
94     ReactionField,
95     Table,
96     TableTwin,
97     Ewald,
98     EwaldTwin,
99     Count
100 };
101
102 //! \brief Lookup function for Coulomb kernel type
103 CoulombKernelType getCoulombKernelType(Nbnxm::EwaldExclusionType ewaldExclusionType,
104                                        CoulombInteractionType    coulombInteractionType,
105                                        bool                      haveEqualCoulombVwdRadii);
106
107 /*! \brief Kinds of Van der Waals treatments in SIMD Verlet kernels
108  *
109  * The \p LJCUT_COMB refers to the LJ combination rule for the short range.
110  * The \p EWALDCOMB refers to the combination rule for the grid part.
111  * \p vdwktNR is the number of VdW treatments for the SIMD kernels.
112  * \p vdwktNR_ref is the number of VdW treatments for the C reference kernels.
113  * These two numbers differ, because currently only the reference kernels
114  * support LB combination rules for the LJ-Ewald grid part.
115  */
116 enum
117 {
118     vdwktLJCUT_COMBGEOM,
119     vdwktLJCUT_COMBLB,
120     vdwktLJCUT_COMBNONE,
121     vdwktLJFORCESWITCH,
122     vdwktLJPOTSWITCH,
123     vdwktLJEWALDCOMBGEOM,
124     vdwktLJEWALDCOMBLB,
125     vdwktNR = vdwktLJEWALDCOMBLB,
126     vdwktNR_ref
127 };
128
129 /*! \brief Clears the force buffer.
130  *
131  * Either the whole buffer is cleared or only the parts used
132  * by thread/task \p outputIndex when nbat->bUseBufferFlags is set.
133  *
134  * \param[in,out] nbat         The Nbnxm atom data
135  * \param[in]     outputIndex  The index of the output object to clear
136  */
137 void clearForceBuffer(nbnxn_atomdata_t* nbat, int outputIndex);
138
139 /*! \brief Clears the shift forces.
140  */
141 void clear_fshift(real* fshift);
142
143 /*! \brief Reduces the collected energy terms over the pair-lists/threads.
144  */
145 void reduce_energies_over_lists(const nbnxn_atomdata_t* nbat, int nlist, real* Vvdw, real* Vc);
146
147 #endif