2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file defines functions for managing threading of listed
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_listed-forces
46 #include "manage-threading.h"
52 #include "gromacs/listed-forces/bonded.h"
53 #include "gromacs/mdlib/forcerec-threading.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/stringutil.h"
58 //! Divides listed interactions over threads
59 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
66 idef->nthreads = nthreads;
68 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
70 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
71 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
74 for (ftype = 0; ftype < F_NRE; ftype++)
76 if (ftype_is_bonded_potential(ftype))
78 nat1 = interaction_function[ftype].nratoms + 1;
80 for (t = 0; t <= nthreads; t++)
82 /* Divide the interactions equally over the threads.
83 * When the different types of bonded interactions
84 * are distributed roughly equally over the threads,
85 * this should lead to well localized output into
86 * the force buffer on each thread.
87 * If this is not the case, a more advanced scheme
88 * (not implemented yet) will do better.
90 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
92 /* Ensure that distance restraint pairs with the same label
93 * end up on the same thread.
94 * This is slighlty tricky code, since the next for iteration
95 * may have an initial il_nr_thread lower than the final value
96 * in the previous iteration, but this will anyhow be increased
97 * to the approriate value again by this while loop.
99 while (ftype == F_DISRES &&
101 il_nr_thread < idef->il[ftype].nr &&
102 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
103 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
105 il_nr_thread += nat1;
108 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
114 //! Construct a reduction mask for which interaction was computed on which thread
116 calc_bonded_reduction_mask(gmx_bitmask_t *mask,
121 int ftype, nb, nat1, nb0, nb1, i, a;
124 for (ftype = 0; ftype < F_NRE; ftype++)
126 if (ftype_is_bonded_potential(ftype))
128 nb = idef->il[ftype].nr;
131 nat1 = interaction_function[ftype].nratoms + 1;
133 /* Divide this interaction equally over the threads.
134 * This is not stored: should match division in calc_bonds.
136 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
137 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
139 for (i = nb0; i < nb1; i += nat1)
141 for (a = 1; a < nat1; a++)
143 bitmask_set_bit(mask, idef->il[ftype].iatoms[i+a]>>shift);
152 /*! \brief We divide the force array in a maximum of BITMASK_SIZE (default 32) blocks.
153 * Minimum force block reduction size is thus 2^6=64.
155 const int maxBlockBits = BITMASK_SIZE;
157 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
162 assert(fr->nthreads >= 1);
164 /* Divide the bonded interaction over the threads */
165 divide_bondeds_over_threads(idef, fr->nthreads);
167 if (fr->nthreads == 1)
175 while (fr->natoms_force > (int)(maxBlockBits*(1U<<fr->red_ashift)))
181 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
185 /* Determine to which blocks each thread's bonded force calculation
186 * contributes. Store this is a mask for each thread.
188 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
189 for (t = 1; t < fr->nthreads; t++)
191 calc_bonded_reduction_mask(&fr->f_t[t].red_mask,
192 idef, fr->red_ashift, t, fr->nthreads);
195 /* Determine the maximum number of blocks we need to reduce over */
198 for (t = 0; t < fr->nthreads; t++)
201 for (b = 0; b < maxBlockBits; b++)
203 if (bitmask_is_set(fr->f_t[t].red_mask, b))
205 fr->red_nblock = std::max(fr->red_nblock, b+1);
211 #if BITMASK_SIZE <= 64 //move into bitmask when it is C++
212 std::string flags = gmx::formatString("%x", fr->f_t[t].red_mask);
214 std::string flags = gmx::formatAndJoin(fr->f_t[t].red_mask,
215 fr->f_t[t].red_mask+BITMASK_ALEN,
216 "", gmx::StringFormatter("%x"));
218 fprintf(debug, "thread %d flags %s count %d\n",
219 t, flags.c_str(), c);
225 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
226 fr->red_nblock, 1<<fr->red_ashift);
227 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
228 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
229 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));