Clean up thread teardown
[gromacs.git] / src / gromacs / listed-forces / listed-internal.h
blob01b11aa255b53c35c5c20d3fc5450544a55682b0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
35 /*! \internal \file
37 * \brief This file contains declarations for functions needed
38 * internally by the module.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_listed-forces
43 #ifndef GMX_LISTED_FORCES_LISTED_INTERNAL_H
44 #define GMX_LISTED_FORCES_LISTED_INTERNAL_H
46 #include "gromacs/math/vectypes.h"
47 #include "gromacs/mdtypes/enerdata.h"
48 #include "gromacs/topology/idef.h"
49 #include "gromacs/topology/ifunc.h"
50 #include "gromacs/utility/bitmask.h"
52 /* We reduce the force array in blocks of 32 atoms. This is large enough
53 * to not cause overhead and 32*sizeof(rvec) is a multiple of the cache-line
54 * size on all systems.
56 static const int reduction_block_size = 32; /**< Force buffer block size in atoms*/
57 static const int reduction_block_bits = 5; /**< log2(reduction_block_size) */
59 /*! \internal \brief struct with output for bonded forces, used per thread */
60 typedef struct
62 rvec4 *f; /**< Force array */
63 int f_nalloc; /**< Allocation size of f */
64 gmx_bitmask_t *mask; /**< Mask for marking which parts of f are filled, working array for constructing mask in bonded_threading_t */
65 int nblock_used; /**< Number of blocks touched by our thread */
66 int *block_index; /**< Index to touched blocks, size nblock_used */
67 int block_nalloc; /**< Allocation size of f (*reduction_block_size), mask_index, mask */
69 rvec *fshift; /**< Shift force array, size SHIFTS */
70 real ener[F_NRE]; /**< Energy array */
71 gmx_grppairener_t grpp; /**< Group pair energy data for pairs */
72 real dvdl[efptNR]; /**< Free-energy dV/dl output */
74 f_thread_t;
76 /*! \internal \brief struct contain all data for bonded force threading */
77 struct bonded_threading_t
79 /* Thread local force and energy data */
80 int nthreads; /**< Number of threads to be used for bondeds */
81 f_thread_t *f_t; /**< Force/energy data per thread, size nthreads */
82 int nblock_used; /**< The number of force blocks to reduce */
83 int *block_index; /**< Index of size nblock_used into mask */
84 gmx_bitmask_t *mask; /**< 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 */
85 int block_nalloc; /**< Allocation size of block_index and mask */
87 bool haveBondeds; /**< true if we have and thus need to reduce bonded forces */
89 /* There are two different ways to distribute the bonded force calculation
90 * over the threads. We dedice which to use based on the number of threads.
92 int max_nthread_uniform; /**< Maximum thread count for uniform distribution of bondeds over threads */
94 int *il_thread_division; /**< Stores the division of work in the t_list over threads.
96 * The division of the normal bonded interactions of threads.
97 * il_thread_division[ftype*(nthreads+1)+t] contains an index
98 * into il[ftype].iatoms; thread th operates on t=th to t=th+1. */
99 int il_thread_division_nalloc; /**< Allocation size of the work division, at least F_NRE*(nthreads+1). */
103 /*! \brief Returns the global topology atom number belonging to local
104 * atom index i.
106 * This function is intended for writing ascii output and returns atom
107 * numbers starting at 1. When global_atom_index=NULL returns i+1.
110 glatnr(const int *global_atom_index, int i);
112 #endif