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,2019, 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"
58 #include "gromacs/listed_forces/gpubonded.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/topology/ifunc.h"
62 #include "gromacs/utility/exceptions.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/gmxassert.h"
65 #include "gromacs/utility/smalloc.h"
67 #include "listed_internal.h"
68 #include "utilities.h"
70 /*! \brief struct for passing all data required for a function type */
72 const t_ilist *il; /**< pointer to t_ilist entry corresponding to ftype */
73 int ftype; /**< the function type index */
74 int nat; /**< nr of atoms involved in a single ftype interaction */
77 /*! \brief Divides listed interactions over threads
79 * This routine attempts to divide all interactions of the numType bondeds
80 * types stored in ild over the threads such that each thread has roughly
81 * equal load and different threads avoid touching the same atoms as much
84 static void divide_bondeds_by_locality(bonded_threading_t *bt,
86 const ilist_data_t *ild)
89 int ind[F_NRE]; /* index into the ild[].il->iatoms */
90 int at_ind[F_NRE]; /* index of the first atom of the interaction at ind */
93 assert(numType <= F_NRE);
96 for (f = 0; f < numType; f++)
98 /* Sum #bondeds*#atoms_per_bond over all bonded types */
99 nat_tot += ild[f].il->nr/(ild[f].nat + 1)*ild[f].nat;
100 /* The start bound for thread 0 is 0 for all interactions */
102 /* Initialize the next atom index array */
103 assert(ild[f].il->nr > 0);
104 at_ind[f] = ild[f].il->iatoms[1];
108 /* Loop over the end bounds of the nthreads threads to determine
109 * which interactions threads 0 to nthreads shall calculate.
111 * NOTE: The cost of these combined loops is #interactions*numType.
112 * This code is running single threaded (difficult to parallelize
113 * over threads). So the relative cost of this function increases
114 * linearly with the number of threads. Since the inner-most loop
115 * is cheap and this is done only at DD repartitioning, the cost should
116 * be negligble. At high thread count many other parts of the code
117 * scale the same way, so it's (currently) not worth improving this.
119 for (t = 1; t <= bt->nthreads; t++)
123 /* Here we assume that the computational cost is proportional
124 * to the number of atoms in the interaction. This is a rough
125 * measure, but roughly correct. Usually there are very few
126 * interactions anyhow and there are distributed relatively
127 * uniformly. Proper and RB dihedrals are often distributed
128 * non-uniformly, but their cost is roughly equal.
130 nat_thread = (nat_tot*t)/bt->nthreads;
132 while (nat_sum < nat_thread)
134 /* To divide bonds based on atom order, we compare
135 * the index of the first atom in the bonded interaction.
136 * This works well, since the domain decomposition generates
137 * bondeds in order of the atoms by looking up interactions
138 * which are linked to the first atom in each interaction.
139 * It usually also works well without DD, since than the atoms
140 * in bonded interactions are usually in increasing order.
141 * If they are not assigned in increasing order, the balancing
142 * is still good, but the memory access and reduction cost will
147 /* Find out which of the types has the lowest atom index */
149 for (f = 1; f < numType; f++)
151 if (at_ind[f] < at_ind[f_min])
156 assert(f_min >= 0 && f_min < numType);
158 /* Assign the interaction with the lowest atom index (of type
159 * index f_min) to thread t-1 by increasing ind.
161 ind[f_min] += ild[f_min].nat + 1;
162 nat_sum += ild[f_min].nat;
164 /* Update the first unassigned atom index for this type */
165 if (ind[f_min] < ild[f_min].il->nr)
167 at_ind[f_min] = ild[f_min].il->iatoms[ind[f_min] + 1];
171 /* We have assigned all interactions of this type.
172 * Setting at_ind to INT_MAX ensures this type will not be
173 * chosen in the for loop above during next iterations.
175 at_ind[f_min] = INT_MAX;
179 /* Store the bonded end boundaries (at index t) for thread t-1 */
180 for (f = 0; f < numType; f++)
182 bt->workDivision.setBound(ild[f].ftype, t, ind[f]);
186 for (f = 0; f < numType; f++)
188 assert(ind[f] == ild[f].il->nr);
192 //! Return whether function type \p ftype in \p idef has perturbed interactions
193 static bool ftypeHasPerturbedEntries(const t_idef &idef,
196 GMX_ASSERT(idef.ilsort == ilsortNO_FE || idef.ilsort == ilsortFE_SORTED,
197 "Perturbed interations should be sorted here");
199 const t_ilist &ilist = idef.il[ftype];
201 return (idef.ilsort != ilsortNO_FE && ilist.nr_nonperturbed != ilist.nr);
204 //! Divides bonded interactions over threads and GPU
205 static void divide_bondeds_over_threads(bonded_threading_t *bt,
206 bool useGpuForBondeds,
209 ilist_data_t ild[F_NRE];
211 assert(bt->nthreads > 0);
213 bt->haveBondeds = false;
215 size_t fTypeGpuIndex = 0;
216 for (int fType = 0; fType < F_NRE; fType++)
218 if (!ftype_is_bonded_potential(fType))
223 const t_ilist &il = idef.il[fType];
224 int nrToAssignToCpuThreads = il.nr;
226 if (useGpuForBondeds &&
227 fTypeGpuIndex < gmx::fTypesOnGpu.size() &&
228 gmx::fTypesOnGpu[fTypeGpuIndex] == fType)
232 /* Perturbation is not implemented in the GPU bonded kernels.
233 * But instead of doing all on the CPU, we could do only
234 * the actually perturbed interactions on the CPU.
236 if (!ftypeHasPerturbedEntries(idef, fType))
238 /* We will assign this interaction type to the GPU */
239 nrToAssignToCpuThreads = 0;
243 if (nrToAssignToCpuThreads > 0)
245 bt->haveBondeds = true;
248 if (nrToAssignToCpuThreads == 0)
250 /* No interactions, avoid all the integer math below */
251 for (int t = 0; t <= bt->nthreads; t++)
253 bt->workDivision.setBound(fType, t, 0);
256 else if (bt->nthreads <= bt->max_nthread_uniform || fType == F_DISRES)
258 /* On up to 4 threads, load balancing the bonded work
259 * is more important than minimizing the reduction cost.
262 const int stride = 1 + NRAL(fType);
264 for (int t = 0; t <= bt->nthreads; t++)
266 /* Divide equally over the threads */
267 int nr_t = (((nrToAssignToCpuThreads/stride)*t)/bt->nthreads)*stride;
269 if (fType == F_DISRES)
271 /* Ensure that distance restraint pairs with the same label
272 * end up on the same thread.
274 while (nr_t > 0 && nr_t < nrToAssignToCpuThreads &&
275 idef.iparams[il.iatoms[nr_t]].disres.label ==
276 idef.iparams[il.iatoms[nr_t - stride]].disres.label)
282 bt->workDivision.setBound(fType, t, nr_t);
287 /* Add this fType to the list to be distributed */
288 int nat = NRAL(fType);
289 ild[numType].ftype = fType;
290 ild[numType].il = &il;
291 ild[numType].nat = nat;
293 /* The first index for the thread division is always 0 */
294 bt->workDivision.setBound(fType, 0, 0);
302 divide_bondeds_by_locality(bt, numType, ild);
309 fprintf(debug, "Division of bondeds over threads:\n");
310 for (f = 0; f < F_NRE; f++)
312 if (ftype_is_bonded_potential(f) && idef.il[f].nr > 0)
316 fprintf(debug, "%16s", interaction_function[f].name);
317 for (t = 0; t < bt->nthreads; t++)
319 fprintf(debug, " %4d",
320 (bt->workDivision.bound(f, t + 1) -
321 bt->workDivision.bound(f, t))/
324 fprintf(debug, "\n");
330 //! Construct a reduction mask for which parts (blocks) of the force array are touched on which thread task
332 calc_bonded_reduction_mask(int natoms,
333 f_thread_t *f_thread,
336 const bonded_threading_t &bondedThreading)
338 static_assert(BITMASK_SIZE == GMX_OPENMP_MAX_THREADS, "For the error message below we assume these two are equal.");
340 if (bondedThreading.nthreads > BITMASK_SIZE)
343 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 passed to CMake.",
344 bondedThreading.nthreads, GMX_OPENMP_MAX_THREADS);
347 GMX_ASSERT(bondedThreading.nthreads <= BITMASK_SIZE, "We need at least nthreads bits in the mask");
349 int nblock = (natoms + reduction_block_size - 1) >> reduction_block_bits;
351 if (nblock > f_thread->block_nalloc)
353 f_thread->block_nalloc = over_alloc_large(nblock);
354 srenew(f_thread->mask, f_thread->block_nalloc);
355 srenew(f_thread->block_index, f_thread->block_nalloc);
356 // NOTE: It seems f_thread->f does not need to be aligned
357 sfree_aligned(f_thread->f);
358 snew_aligned(f_thread->f, f_thread->block_nalloc*reduction_block_size, 128);
361 gmx_bitmask_t *mask = f_thread->mask;
363 for (int b = 0; b < nblock; b++)
365 bitmask_clear(&mask[b]);
368 for (int ftype = 0; ftype < F_NRE; ftype++)
370 if (ftype_is_bonded_potential(ftype))
372 int nb = idef.il[ftype].nr;
375 int nat1 = interaction_function[ftype].nratoms + 1;
377 int nb0 = bondedThreading.workDivision.bound(ftype, thread);
378 int nb1 = bondedThreading.workDivision.bound(ftype, thread + 1);
380 for (int i = nb0; i < nb1; i += nat1)
382 for (int a = 1; a < nat1; a++)
384 bitmask_set_bit(&mask[idef.il[ftype].iatoms[i+a] >> reduction_block_bits], thread);
391 /* Make an index of the blocks our thread touches, so we can do fast
392 * force buffer clearing.
394 f_thread->nblock_used = 0;
395 for (int b = 0; b < nblock; b++)
397 if (bitmask_is_set(mask[b], thread))
399 f_thread->block_index[f_thread->nblock_used++] = b;
404 void setup_bonded_threading(bonded_threading_t *bt,
406 bool useGpuForBondeds,
411 assert(bt->nthreads >= 1);
413 /* Divide the bonded interaction over the threads */
414 divide_bondeds_over_threads(bt, useGpuForBondeds, idef);
416 if (!bt->haveBondeds)
418 /* We don't have bondeds, so there is nothing to reduce */
422 /* Determine to which blocks each thread's bonded force calculation
423 * contributes. Store this as a mask for each thread.
425 #pragma omp parallel for num_threads(bt->nthreads) schedule(static)
426 for (int t = 0; t < bt->nthreads; t++)
430 calc_bonded_reduction_mask(numAtoms, bt->f_t[t].get(),
433 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
436 /* Reduce the masks over the threads and determine which blocks
437 * we need to reduce over.
439 int nblock_tot = (numAtoms + reduction_block_size - 1) >> reduction_block_bits;
440 /* Ensure we have sufficient space for all blocks */
441 if (static_cast<size_t>(nblock_tot) > bt->block_index.size())
443 bt->block_index.resize(nblock_tot);
445 if (static_cast<size_t>(nblock_tot) > bt->mask.size())
447 bt->mask.resize(nblock_tot);
450 for (int b = 0; b < nblock_tot; b++)
452 gmx_bitmask_t *mask = &bt->mask[b];
454 /* Generate the union over the threads of the bitmask */
456 for (int t = 0; t < bt->nthreads; t++)
458 bitmask_union(mask, bt->f_t[t]->mask[b]);
460 if (!bitmask_is_zero(*mask))
462 bt->block_index[bt->nblock_used++] = b;
468 for (int t = 0; t < bt->nthreads; t++)
470 if (bitmask_is_set(*mask, t))
479 fprintf(debug, "block %d flags %s count %d\n",
480 b, to_hex_string(*mask).c_str(), c);
486 fprintf(debug, "Number of %d atom blocks to reduce: %d\n",
487 reduction_block_size, bt->nblock_used);
488 fprintf(debug, "Reduction density %.2f for touched blocks only %.2f\n",
489 ctot*reduction_block_size/static_cast<double>(numAtoms),
490 ctot/static_cast<double>(bt->nblock_used));
494 void tear_down_bonded_threading(bonded_threading_t *bt)
499 f_thread_t::f_thread_t(int numEnergyGroups) :
500 grpp(numEnergyGroups)
502 snew(fshift, SHIFTS);
505 f_thread_t::~f_thread_t()
513 bonded_threading_t::bonded_threading_t(const int numThreads,
514 const int numEnergyGroups) :
515 nthreads(numThreads),
518 workDivision(nthreads),
519 foreignLambdaWorkDivision(1)
521 f_t.resize(numThreads);
522 #pragma omp parallel for num_threads(nthreads) schedule(static)
523 for (int t = 0; t < nthreads; t++)
527 /* Note that thread 0 uses the global fshift and energy arrays,
528 * but to keep the code simple, we initialize all data here.
530 f_t[t] = std::make_unique<f_thread_t>(numEnergyGroups);
532 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
536 bonded_threading_t *init_bonded_threading(FILE *fplog,
539 /* These thread local data structures are used for bondeds only.
541 * Note that we also use there structures when running single-threaded.
542 * This is because the bonded force buffer uses type rvec4, whereas
543 * the normal force buffer is uses type rvec. This leads to a little
544 * reduction overhead, but the speed gain in the bonded calculations
545 * of doing transposeScatterIncr/DecrU with aligment 4 instead of 3
546 * is much larger than the reduction overhead.
548 bonded_threading_t *bt = new bonded_threading_t(gmx_omp_nthreads_get(emntBonded),
551 /* The optimal value after which to switch from uniform to localized
552 * bonded interaction distribution is 3, 4 or 5 depending on the system
555 const int max_nthread_uniform = 4;
558 if ((ptr = getenv("GMX_BONDED_NTHREAD_UNIFORM")) != nullptr)
560 sscanf(ptr, "%d", &bt->max_nthread_uniform);
561 if (fplog != nullptr)
563 fprintf(fplog, "\nMax threads for uniform bonded distribution set to %d by env.var.\n",
564 bt->max_nthread_uniform);
569 bt->max_nthread_uniform = max_nthread_uniform;