Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / listed_forces / listed_internal.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,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.
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 /*! \internal \file
36  *
37  * \brief This file contains declarations for functions needed
38  * internally by the module.
39  *
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_listed_forces
42  */
43 #ifndef GMX_LISTED_FORCES_LISTED_INTERNAL_H
44 #define GMX_LISTED_FORCES_LISTED_INTERNAL_H
45
46 #include <memory>
47
48 #include "gromacs/math/vectypes.h"
49 #include "gromacs/mdtypes/enerdata.h"
50 #include "gromacs/topology/idef.h"
51 #include "gromacs/topology/ifunc.h"
52 #include "gromacs/utility/bitmask.h"
53
54 /* We reduce the force array in blocks of 32 atoms. This is large enough
55  * to not cause overhead and 32*sizeof(rvec) is a multiple of the cache-line
56  * size on all systems.
57  */
58 static const int reduction_block_size = 32; /**< Force buffer block size in atoms*/
59 static const int reduction_block_bits = 5;  /**< log2(reduction_block_size) */
60
61 /*! \internal \brief The division of bonded interactions of the threads */
62 class WorkDivision
63 {
64 public:
65     //! Constructor
66     WorkDivision(int numThreads) : stride_(numThreads + 1), packedBounds_(F_NRE * stride_) {}
67
68     //! Sets the bound between threads \p boundIndex-1 and \p boundIndex to \p count
69     void setBound(int functionType, int boundIndex, int count)
70     {
71         packedBounds_[functionType * stride_ + boundIndex] = count;
72     }
73
74     //! Returns the bound between threads \p boundIndex-1 and \p boundIndex
75     int bound(int functionType, int boundIndex) const
76     {
77         return packedBounds_[functionType * stride_ + boundIndex];
78     }
79
80     //! Returns the last bound
81     int end(int ftype) const { return bound(ftype, stride_ - 1); }
82
83 private:
84     //! The stride_ between and size of the entries for a function type
85     int stride_;
86     //! The bounds stored as a flat array for fast access
87     std::vector<int> packedBounds_;
88 };
89
90 /*! \internal \brief struct with output for bonded forces, used per thread */
91 struct f_thread_t
92 {
93     //! Constructor
94     f_thread_t(int numEnergyGroups);
95
96     ~f_thread_t();
97
98     rvec4*         f        = nullptr; /**< Force array */
99     int            f_nalloc = 0;       /**< Allocation size of f */
100     gmx_bitmask_t* mask =
101             nullptr; /**< Mask for marking which parts of f are filled, working array for constructing mask in bonded_threading_t */
102     int  nblock_used;            /**< Number of blocks touched by our thread */
103     int* block_index  = nullptr; /**< Index to touched blocks, size nblock_used */
104     int  block_nalloc = 0; /**< Allocation size of f (*reduction_block_size), mask_index, mask */
105
106     rvec*             fshift;       /**< Shift force array, size SHIFTS */
107     real              ener[F_NRE];  /**< Energy array */
108     gmx_grppairener_t grpp;         /**< Group pair energy data for pairs */
109     real              dvdl[efptNR]; /**< Free-energy dV/dl output */
110 };
111
112 /*! \internal \brief struct contain all data for bonded force threading */
113 struct bonded_threading_t
114 {
115     //! Constructor
116     bonded_threading_t(int numThreads, int numEnergyGroups);
117
118     //! Number of threads to be used for bondeds
119     int nthreads;
120     //! Force/energy data per thread, size nthreads, stored in unique_ptr to allow thread local allocation
121     std::vector<std::unique_ptr<f_thread_t>> f_t;
122     //! The number of force blocks to reduce
123     int nblock_used;
124     //! Index of size nblock_used into mask
125     std::vector<int> block_index;
126     //! Mask array, one element corresponds to a block of reduction_block_size atoms of the force array, bit corresponding to thread indices set if a thread writes to that block
127     std::vector<gmx_bitmask_t> mask;
128     //! true if we have and thus need to reduce bonded forces
129     bool haveBondeds;
130
131     /* There are two different ways to distribute the bonded force calculation
132      * over the threads. We dedice which to use based on the number of threads.
133      */
134     //! Maximum thread count for uniform distribution of bondeds over threads
135     int max_nthread_uniform;
136
137     //! The division of work in the t_list over threads.
138     WorkDivision workDivision;
139
140     //! Work division for free-energy foreign lambda calculations, always uses 1 thread
141     WorkDivision foreignLambdaWorkDivision;
142 };
143
144
145 /*! \brief Returns the global topology atom number belonging to local
146  * atom index i.
147  *
148  * This function is intended for writing ascii output and returns atom
149  * numbers starting at 1.  When global_atom_index=NULL returns i+1.
150  */
151 int glatnr(const int* global_atom_index, int i);
152
153 #endif