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,2015, 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"
53 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
54 #include "gromacs/listed-forces/listed-forces.h"
55 #include "gromacs/mdlib/forcerec-threading.h"
56 #include "gromacs/pbcutil/ishift.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
59 #include "gromacs/utility/stringutil.h"
61 /*! \brief struct for passing all data required for a function type */
63 int ftype; /**< the function type index */
64 t_ilist *il; /**< pointer to t_ilist entry corresponding to ftype */
65 int nat; /**< nr of atoms involved in a single ftype interaction */
68 /*! \brief Divides listed interactions over threads
70 * This routine attempts to divide all interactions of the ntype bondeds
71 * types stored in ild over the threads such that each thread has roughly
72 * equal load and different threads avoid touching the same atoms as much
75 static void divide_bondeds_by_locality(int ntype,
76 const ilist_data_t *ild,
81 int ind[F_NRE]; /* index into the ild[].il->iatoms */
82 int at_ind[F_NRE]; /* index of the first atom of the interaction at ind */
85 assert(ntype <= F_NRE);
88 for (f = 0; f < ntype; f++)
90 /* Sum #bondeds*#atoms_per_bond over all bonded types */
91 nat_tot += ild[f].il->nr/(ild[f].nat + 1)*ild[f].nat;
92 /* The start bound for thread 0 is 0 for all interactions */
94 /* Initialize the next atom index array */
95 assert(ild[f].il->nr > 0);
96 at_ind[f] = ild[f].il->iatoms[1];
100 /* Loop over the end bounds of the nthread threads to determine
101 * which interactions threads 0 to nthread shall calculate.
103 * NOTE: The cost of these combined loops is #interactions*ntype.
104 * This code is running single threaded (difficult to parallelize
105 * over threads). So the relative cost of this function increases
106 * linearly with the number of threads. Since the inner-most loop
107 * is cheap and this is done only at DD repartitioning, the cost should
108 * be negligble. At high thread count many other parts of the code
109 * scale the same way, so it's (currently) not worth improving this.
111 for (t = 1; t <= nthread; t++)
115 /* Here we assume that the computational cost is proportional
116 * to the number of atoms in the interaction. This is a rough
117 * measure, but roughly correct. Usually there are very few
118 * interactions anyhow and there are distributed relatively
119 * uniformly. Proper and RB dihedrals are often distributed
120 * non-uniformly, but their cost is roughly equal.
122 nat_thread = (nat_tot*t)/nthread;
124 while (nat_sum < nat_thread)
126 /* To divide bonds based on atom order, we compare
127 * the index of the first atom in the bonded interaction.
128 * This works well, since the domain decomposition generates
129 * bondeds in order of the atoms by looking up interactions
130 * which are linked to the first atom in each interaction.
131 * It usually also works well without DD, since than the atoms
132 * in bonded interactions are usually in increasing order.
133 * If they are not assigned in increasing order, the balancing
134 * is still good, but the memory access and reduction cost will
139 /* Find out which of the types has the lowest atom index */
141 for (f = 1; f < ntype; f++)
143 if (at_ind[f] < at_ind[f_min])
148 assert(f_min >= 0 && f_min < ntype);
150 /* Assign the interaction with the lowest atom index (of type
151 * index f_min) to thread t-1 by increasing ind.
153 ind[f_min] += ild[f_min].nat + 1;
154 nat_sum += ild[f_min].nat;
156 /* Update the first unassigned atom index for this type */
157 if (ind[f_min] < ild[f_min].il->nr)
159 at_ind[f_min] = ild[f_min].il->iatoms[ind[f_min] + 1];
163 /* We have assigned all interactions of this type.
164 * Setting at_ind to INT_MAX ensures this type will not be
165 * chosen in the for loop above during next iterations.
167 at_ind[f_min] = INT_MAX;
171 /* Store the bonded end boundaries (at index t) for thread t-1 */
172 for (f = 0; f < ntype; f++)
174 idef->il_thread_division[ild[f].ftype*(nthread + 1) + t] = ind[f];
178 for (f = 0; f < ntype; f++)
180 assert(ind[f] == ild[f].il->nr);
184 //! Divides bonded interactions over threads
185 static void divide_bondeds_over_threads(t_idef *idef,
187 int max_nthread_uniform)
189 ilist_data_t ild[F_NRE];
195 idef->nthreads = nthread;
197 if (F_NRE*(nthread + 1) > idef->il_thread_division_nalloc)
199 idef->il_thread_division_nalloc = F_NRE*(nthread + 1);
200 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
204 for (f = 0; f < F_NRE; f++)
206 if (!ftype_is_bonded_potential(f))
211 if (idef->il[f].nr == 0)
213 /* No interactions, avoid all the integer math below */
215 for (t = 0; t <= nthread; t++)
217 idef->il_thread_division[f*(nthread + 1) + t] = 0;
220 else if (nthread <= max_nthread_uniform || f == F_DISRES)
222 /* On up to 4 threads, load balancing the bonded work
223 * is more important than minimizing the reduction cost.
227 /* nat1 = 1 + #atoms(ftype) which is the stride use for iatoms */
230 for (t = 0; t <= nthread; t++)
234 /* Divide equally over the threads */
235 nr_t = (((idef->il[f].nr/nat1)*t)/nthread)*nat1;
239 /* Ensure that distance restraint pairs with the same label
240 * end up on the same thread.
242 while (nr_t > 0 && nr_t < idef->il[f].nr &&
243 idef->iparams[idef->il[f].iatoms[nr_t]].disres.label ==
244 idef->iparams[idef->il[f].iatoms[nr_t-nat1]].disres.label)
250 idef->il_thread_division[f*(nthread + 1) + t] = nr_t;
255 /* Add this ftype to the list to be distributed */
259 ild[ntype].ftype = f;
260 ild[ntype].il = &idef->il[f];
261 ild[ntype].nat = nat;
263 /* The first index for the thread division is always 0 */
264 idef->il_thread_division[f*(nthread + 1)] = 0;
272 divide_bondeds_by_locality(ntype, ild, nthread, idef);
279 fprintf(debug, "Division of bondeds over threads:\n");
280 for (f = 0; f < F_NRE; f++)
282 if (ftype_is_bonded_potential(f) && idef->il[f].nr > 0)
286 fprintf(debug, "%16s", interaction_function[f].name);
287 for (t = 0; t < nthread; t++)
289 fprintf(debug, " %4d",
290 (idef->il_thread_division[f*(nthread + 1) + t + 1] -
291 idef->il_thread_division[f*(nthread + 1) + t])/
294 fprintf(debug, "\n");
300 //! Construct a reduction mask for which interaction was computed on which thread
302 calc_bonded_reduction_mask(gmx_bitmask_t *mask,
307 int ftype, nb, nat1, nb0, nb1, i, a;
310 for (ftype = 0; ftype < F_NRE; ftype++)
312 if (ftype_is_bonded_potential(ftype))
314 nb = idef->il[ftype].nr;
317 nat1 = interaction_function[ftype].nratoms + 1;
319 /* Divide this interaction equally over the threads.
320 * This is not stored: should match division in calc_bonds.
322 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
323 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
325 for (i = nb0; i < nb1; i += nat1)
327 for (a = 1; a < nat1; a++)
329 bitmask_set_bit(mask, idef->il[ftype].iatoms[i+a]>>shift);
338 /*! \brief We divide the force array in a maximum of BITMASK_SIZE (default 32) blocks.
339 * Minimum force block reduction size is thus 2^6=64.
341 const int maxBlockBits = BITMASK_SIZE;
343 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
348 assert(fr->nthreads >= 1);
350 /* Divide the bonded interaction over the threads */
351 divide_bondeds_over_threads(idef,
353 fr->bonded_max_nthread_uniform);
355 if (fr->nthreads == 1)
363 while (fr->natoms_force > (int)(maxBlockBits*(1U<<fr->red_ashift)))
369 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
373 /* Determine to which blocks each thread's bonded force calculation
374 * contributes. Store this is a mask for each thread.
376 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
377 for (t = 1; t < fr->nthreads; t++)
379 calc_bonded_reduction_mask(&fr->f_t[t].red_mask,
380 idef, fr->red_ashift, t, fr->nthreads);
383 /* Determine the maximum number of blocks we need to reduce over */
386 for (t = 0; t < fr->nthreads; t++)
389 for (b = 0; b < maxBlockBits; b++)
391 if (bitmask_is_set(fr->f_t[t].red_mask, b))
393 fr->red_nblock = std::max(fr->red_nblock, b+1);
399 #if BITMASK_SIZE <= 64 //move into bitmask when it is C++
400 std::string flags = gmx::formatString("%x", fr->f_t[t].red_mask);
402 std::string flags = gmx::formatAndJoin(fr->f_t[t].red_mask,
403 fr->f_t[t].red_mask+BITMASK_ALEN,
404 "", gmx::StringFormatter("%x"));
406 fprintf(debug, "thread %d flags %s count %d\n",
407 t, flags.c_str(), c);
413 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
414 fr->red_nblock, 1<<fr->red_ashift);
415 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
416 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
417 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
421 void init_bonded_threading(FILE *fplog, t_forcerec *fr, int nenergrp)
423 /* These thread local data structures are used for bondeds only */
424 fr->nthreads = gmx_omp_nthreads_get(emntBonded);
426 if (fr->nthreads > 1)
430 snew(fr->f_t, fr->nthreads);
431 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
432 for (t = 0; t < fr->nthreads; t++)
434 /* Thread 0 uses the global force and energy arrays */
440 fr->f_t[t].f_nalloc = 0;
441 snew(fr->f_t[t].fshift, SHIFTS);
442 fr->f_t[t].grpp.nener = nenergrp*nenergrp;
443 for (i = 0; i < egNR; i++)
445 snew(fr->f_t[t].grpp.ener[i], fr->f_t[t].grpp.nener);
450 /* The optimal value after which to switch from uniform to localized
451 * bonded interaction distribution is 3, 4 or 5 depending on the system
454 const int max_nthread_uniform = 4;
457 if ((ptr = getenv("GMX_BONDED_NTHREAD_UNIFORM")) != NULL)
459 sscanf(ptr, "%d", &fr->bonded_max_nthread_uniform);
462 fprintf(fplog, "\nMax threads for uniform bonded distribution set to %d by env.var.\n",
463 fr->bonded_max_nthread_uniform);
468 fr->bonded_max_nthread_uniform = max_nthread_uniform;