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,2016,2017,2018, 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"
56 #include "gromacs/listed-forces/listed-forces.h"
57 #include "gromacs/mdlib/gmx_omp_nthreads.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/topology/ifunc.h"
60 #include "gromacs/utility/exceptions.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "gromacs/utility/stringutil.h"
65 #include "listed-internal.h"
67 /*! \brief struct for passing all data required for a function type */
69 const t_ilist *il; /**< pointer to t_ilist entry corresponding to ftype */
70 int ftype; /**< the function type index */
71 int nat; /**< nr of atoms involved in a single ftype interaction */
74 /*! \brief Divides listed interactions over threads
76 * This routine attempts to divide all interactions of the ntype bondeds
77 * types stored in ild over the threads such that each thread has roughly
78 * equal load and different threads avoid touching the same atoms as much
81 static void divide_bondeds_by_locality(bonded_threading_t *bt,
83 const ilist_data_t *ild)
86 int ind[F_NRE]; /* index into the ild[].il->iatoms */
87 int at_ind[F_NRE]; /* index of the first atom of the interaction at ind */
90 assert(ntype <= F_NRE);
93 for (f = 0; f < ntype; f++)
95 /* Sum #bondeds*#atoms_per_bond over all bonded types */
96 nat_tot += ild[f].il->nr/(ild[f].nat + 1)*ild[f].nat;
97 /* The start bound for thread 0 is 0 for all interactions */
99 /* Initialize the next atom index array */
100 assert(ild[f].il->nr > 0);
101 at_ind[f] = ild[f].il->iatoms[1];
105 /* Loop over the end bounds of the nthreads threads to determine
106 * which interactions threads 0 to nthreads shall calculate.
108 * NOTE: The cost of these combined loops is #interactions*ntype.
109 * This code is running single threaded (difficult to parallelize
110 * over threads). So the relative cost of this function increases
111 * linearly with the number of threads. Since the inner-most loop
112 * is cheap and this is done only at DD repartitioning, the cost should
113 * be negligble. At high thread count many other parts of the code
114 * scale the same way, so it's (currently) not worth improving this.
116 for (t = 1; t <= bt->nthreads; t++)
120 /* Here we assume that the computational cost is proportional
121 * to the number of atoms in the interaction. This is a rough
122 * measure, but roughly correct. Usually there are very few
123 * interactions anyhow and there are distributed relatively
124 * uniformly. Proper and RB dihedrals are often distributed
125 * non-uniformly, but their cost is roughly equal.
127 nat_thread = (nat_tot*t)/bt->nthreads;
129 while (nat_sum < nat_thread)
131 /* To divide bonds based on atom order, we compare
132 * the index of the first atom in the bonded interaction.
133 * This works well, since the domain decomposition generates
134 * bondeds in order of the atoms by looking up interactions
135 * which are linked to the first atom in each interaction.
136 * It usually also works well without DD, since than the atoms
137 * in bonded interactions are usually in increasing order.
138 * If they are not assigned in increasing order, the balancing
139 * is still good, but the memory access and reduction cost will
144 /* Find out which of the types has the lowest atom index */
146 for (f = 1; f < ntype; f++)
148 if (at_ind[f] < at_ind[f_min])
153 assert(f_min >= 0 && f_min < ntype);
155 /* Assign the interaction with the lowest atom index (of type
156 * index f_min) to thread t-1 by increasing ind.
158 ind[f_min] += ild[f_min].nat + 1;
159 nat_sum += ild[f_min].nat;
161 /* Update the first unassigned atom index for this type */
162 if (ind[f_min] < ild[f_min].il->nr)
164 at_ind[f_min] = ild[f_min].il->iatoms[ind[f_min] + 1];
168 /* We have assigned all interactions of this type.
169 * Setting at_ind to INT_MAX ensures this type will not be
170 * chosen in the for loop above during next iterations.
172 at_ind[f_min] = INT_MAX;
176 /* Store the bonded end boundaries (at index t) for thread t-1 */
177 for (f = 0; f < ntype; f++)
179 bt->il_thread_division[ild[f].ftype*(bt->nthreads + 1) + t] = ind[f];
183 for (f = 0; f < ntype; f++)
185 assert(ind[f] == ild[f].il->nr);
189 //! Divides bonded interactions over threads
190 static void divide_bondeds_over_threads(bonded_threading_t *bt,
193 ilist_data_t ild[F_NRE];
197 assert(bt->nthreads > 0);
199 if (F_NRE*(bt->nthreads + 1) > bt->il_thread_division_nalloc)
201 bt->il_thread_division_nalloc = F_NRE*(bt->nthreads + 1);
202 srenew(bt->il_thread_division, bt->il_thread_division_nalloc);
205 bt->haveBondeds = false;
207 for (f = 0; f < F_NRE; f++)
209 if (!ftype_is_bonded_potential(f))
214 if (idef->il[f].nr > 0)
216 bt->haveBondeds = true;
219 if (idef->il[f].nr == 0)
221 /* No interactions, avoid all the integer math below */
223 for (t = 0; t <= bt->nthreads; t++)
225 bt->il_thread_division[f*(bt->nthreads + 1) + t] = 0;
228 else if (bt->nthreads <= bt->max_nthread_uniform || f == F_DISRES)
230 /* On up to 4 threads, load balancing the bonded work
231 * is more important than minimizing the reduction cost.
235 /* nat1 = 1 + #atoms(ftype) which is the stride use for iatoms */
238 for (t = 0; t <= bt->nthreads; t++)
242 /* Divide equally over the threads */
243 nr_t = (((idef->il[f].nr/nat1)*t)/bt->nthreads)*nat1;
247 /* Ensure that distance restraint pairs with the same label
248 * end up on the same thread.
250 while (nr_t > 0 && nr_t < idef->il[f].nr &&
251 idef->iparams[idef->il[f].iatoms[nr_t]].disres.label ==
252 idef->iparams[idef->il[f].iatoms[nr_t-nat1]].disres.label)
258 bt->il_thread_division[f*(bt->nthreads + 1) + t] = nr_t;
263 /* Add this ftype to the list to be distributed */
267 ild[ntype].ftype = f;
268 ild[ntype].il = &idef->il[f];
269 ild[ntype].nat = nat;
271 /* The first index for the thread division is always 0 */
272 bt->il_thread_division[f*(bt->nthreads + 1)] = 0;
280 divide_bondeds_by_locality(bt, ntype, ild);
287 fprintf(debug, "Division of bondeds over threads:\n");
288 for (f = 0; f < F_NRE; f++)
290 if (ftype_is_bonded_potential(f) && idef->il[f].nr > 0)
294 fprintf(debug, "%16s", interaction_function[f].name);
295 for (t = 0; t < bt->nthreads; t++)
297 fprintf(debug, " %4d",
298 (bt->il_thread_division[f*(bt->nthreads + 1) + t + 1] -
299 bt->il_thread_division[f*(bt->nthreads + 1) + t])/
302 fprintf(debug, "\n");
308 //! Construct a reduction mask for which parts (blocks) of the force array are touched on which thread task
310 calc_bonded_reduction_mask(int natoms,
311 f_thread_t *f_thread,
314 const bonded_threading_t &bondedThreading)
316 static_assert(BITMASK_SIZE == GMX_OPENMP_MAX_THREADS, "For the error message below we assume these two are equal.");
318 if (bondedThreading.nthreads > BITMASK_SIZE)
321 gmx_fatal(FARGS, "You are using %d OpenMP threads, which is larger than GMX_OPENMP_MAX_THREADS (%d). Decrease the number of OpenMP threads or rebuild GROMACS with a larger value for GMX_OPENMP_MAX_THREADS.",
322 bondedThreading.nthreads, GMX_OPENMP_MAX_THREADS);
325 GMX_ASSERT(bondedThreading.nthreads <= BITMASK_SIZE, "We need at least nthreads bits in the mask");
327 int nblock = (natoms + reduction_block_size - 1) >> reduction_block_bits;
329 if (nblock > f_thread->block_nalloc)
331 f_thread->block_nalloc = over_alloc_large(nblock);
332 srenew(f_thread->mask, f_thread->block_nalloc);
333 srenew(f_thread->block_index, f_thread->block_nalloc);
334 sfree_aligned(f_thread->f);
335 snew_aligned(f_thread->f, f_thread->block_nalloc*reduction_block_size, 128);
338 gmx_bitmask_t *mask = f_thread->mask;
340 for (int b = 0; b < nblock; b++)
342 bitmask_clear(&mask[b]);
345 for (int ftype = 0; ftype < F_NRE; ftype++)
347 if (ftype_is_bonded_potential(ftype))
349 int nb = idef->il[ftype].nr;
352 int nat1 = interaction_function[ftype].nratoms + 1;
354 int nb0 = bondedThreading.il_thread_division[ftype*(bondedThreading.nthreads + 1) + thread];
355 int nb1 = bondedThreading.il_thread_division[ftype*(bondedThreading.nthreads + 1) + thread + 1];
357 for (int i = nb0; i < nb1; i += nat1)
359 for (int a = 1; a < nat1; a++)
361 bitmask_set_bit(&mask[idef->il[ftype].iatoms[i+a] >> reduction_block_bits], thread);
368 /* Make an index of the blocks our thread touches, so we can do fast
369 * force buffer clearing.
371 f_thread->nblock_used = 0;
372 for (int b = 0; b < nblock; b++)
374 if (bitmask_is_set(mask[b], thread))
376 f_thread->block_index[f_thread->nblock_used++] = b;
381 void setup_bonded_threading(bonded_threading_t *bt,
387 assert(bt->nthreads >= 1);
389 /* Divide the bonded interaction over the threads */
390 divide_bondeds_over_threads(bt, idef);
392 if (!bt->haveBondeds)
394 /* We don't have bondeds, so there is nothing to reduce */
398 /* Determine to which blocks each thread's bonded force calculation
399 * contributes. Store this as a mask for each thread.
401 #pragma omp parallel for num_threads(bt->nthreads) schedule(static)
402 for (int t = 0; t < bt->nthreads; t++)
406 calc_bonded_reduction_mask(numAtoms, &bt->f_t[t],
409 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
412 /* Reduce the masks over the threads and determine which blocks
413 * we need to reduce over.
415 int nblock_tot = (numAtoms + reduction_block_size - 1) >> reduction_block_bits;
416 if (nblock_tot > bt->block_nalloc)
418 bt->block_nalloc = over_alloc_large(nblock_tot);
419 srenew(bt->block_index, bt->block_nalloc);
420 srenew(bt->mask, bt->block_nalloc);
423 for (int b = 0; b < nblock_tot; b++)
425 gmx_bitmask_t *mask = &bt->mask[b];
427 /* Generate the union over the threads of the bitmask */
429 for (int t = 0; t < bt->nthreads; t++)
431 bitmask_union(mask, bt->f_t[t].mask[b]);
433 if (!bitmask_is_zero(*mask))
435 bt->block_index[bt->nblock_used++] = b;
441 for (int t = 0; t < bt->nthreads; t++)
443 if (bitmask_is_set(*mask, t))
452 #if BITMASK_SIZE <= 64 //move into bitmask when it is C++
453 std::string flags = gmx::formatString("%x", *mask);
455 std::string flags = gmx::formatAndJoin(*mask,
457 "", gmx::StringFormatter("%x"));
459 fprintf(debug, "block %d flags %s count %d\n",
460 b, flags.c_str(), c);
466 fprintf(debug, "Number of %d atom blocks to reduce: %d\n",
467 reduction_block_size, bt->nblock_used);
468 fprintf(debug, "Reduction density %.2f for touched blocks only %.2f\n",
469 ctot*reduction_block_size/(double)numAtoms,
470 ctot/(double)bt->nblock_used);
474 void tear_down_bonded_threading(bonded_threading_t *bt)
476 for (int th = 0; th < bt->nthreads; th++)
478 sfree(bt->f_t[th].fshift);
479 for (int i = 0; i < egNR; i++)
481 sfree(bt->f_t[th].grpp.ener[i]);
485 sfree(bt->il_thread_division);
489 void init_bonded_threading(FILE *fplog, int nenergrp,
490 struct bonded_threading_t **bt_ptr)
492 bonded_threading_t *bt;
496 /* These thread local data structures are used for bondeds only.
498 * Note that we also use there structures when running single-threaded.
499 * This is because the bonded force buffer uses type rvec4, whereas
500 * the normal force buffer is uses type rvec. This leads to a little
501 * reduction overhead, but the speed gain in the bonded calculations
502 * of doing transposeScatterIncr/DecrU with aligment 4 instead of 3
503 * is much larger than the reduction overhead.
505 bt->nthreads = gmx_omp_nthreads_get(emntBonded);
507 snew(bt->f_t, bt->nthreads);
508 #pragma omp parallel for num_threads(bt->nthreads) schedule(static)
509 for (int t = 0; t < bt->nthreads; t++)
513 /* Note that thread 0 uses the global fshift and energy arrays,
514 * but to keep the code simple, we initialize all data here.
516 bt->f_t[t].f = nullptr;
517 bt->f_t[t].f_nalloc = 0;
518 snew(bt->f_t[t].fshift, SHIFTS);
519 bt->f_t[t].grpp.nener = nenergrp*nenergrp;
520 for (int i = 0; i < egNR; i++)
522 snew(bt->f_t[t].grpp.ener[i], bt->f_t[t].grpp.nener);
525 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
529 bt->block_index = nullptr;
531 bt->block_nalloc = 0;
533 /* The optimal value after which to switch from uniform to localized
534 * bonded interaction distribution is 3, 4 or 5 depending on the system
537 const int max_nthread_uniform = 4;
540 if ((ptr = getenv("GMX_BONDED_NTHREAD_UNIFORM")) != nullptr)
542 sscanf(ptr, "%d", &bt->max_nthread_uniform);
543 if (fplog != nullptr)
545 fprintf(fplog, "\nMax threads for uniform bonded distribution set to %d by env.var.\n",
546 bt->max_nthread_uniform);
551 bt->max_nthread_uniform = max_nthread_uniform;