Move responsibility for bonded threading decomposition
[alexxy/gromacs.git] / src / gromacs / listed-forces / manage-threading.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \internal \file
38  * \brief This file defines functions for managing threading of listed
39  * interactions.
40  *
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_listed-forces
43  */
44 #include "gmxpre.h"
45
46 #include "manage-threading.h"
47
48 #include "config.h"
49
50 #include <assert.h>
51 #include <limits.h>
52 #include <stdlib.h>
53
54 #include <algorithm>
55
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"
64
65 #include "listed-internal.h"
66
67 /*! \brief struct for passing all data required for a function type */
68 typedef struct {
69     int            ftype; /**< the function type index */
70     const t_ilist *il;    /**< pointer to t_ilist entry corresponding to ftype */
71     int            nat;   /**< nr of atoms involved in a single ftype interaction */
72 } ilist_data_t;
73
74 /*! \brief Divides listed interactions over threads
75  *
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
79  * as possible.
80  */
81 static void divide_bondeds_by_locality(bonded_threading_t *bt,
82                                        int                 ntype,
83                                        const ilist_data_t *ild)
84 {
85     int nat_tot, nat_sum;
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 */
88     int f, t;
89
90     assert(ntype <= F_NRE);
91
92     nat_tot = 0;
93     for (f = 0; f < ntype; f++)
94     {
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 */
98         ind[f]    = 0;
99         /* Initialize the next atom index array */
100         assert(ild[f].il->nr > 0);
101         at_ind[f] = ild[f].il->iatoms[1];
102     }
103
104     nat_sum = 0;
105     /* Loop over the end bounds of the nthreads threads to determine
106      * which interactions threads 0 to nthreads shall calculate.
107      *
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.
115      */
116     for (t = 1; t <= bt->nthreads; t++)
117     {
118         int nat_thread;
119
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.
126          */
127         nat_thread = (nat_tot*t)/bt->nthreads;
128
129         while (nat_sum < nat_thread)
130         {
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
140              * be higher.
141              */
142             int f_min;
143
144             /* Find out which of the types has the lowest atom index */
145             f_min = 0;
146             for (f = 1; f < ntype; f++)
147             {
148                 if (at_ind[f] < at_ind[f_min])
149                 {
150                     f_min = f;
151                 }
152             }
153             assert(f_min >= 0 && f_min < ntype);
154
155             /* Assign the interaction with the lowest atom index (of type
156              * index f_min) to thread t-1 by increasing ind.
157              */
158             ind[f_min] += ild[f_min].nat + 1;
159             nat_sum    += ild[f_min].nat;
160
161             /* Update the first unassigned atom index for this type */
162             if (ind[f_min] < ild[f_min].il->nr)
163             {
164                 at_ind[f_min] = ild[f_min].il->iatoms[ind[f_min] + 1];
165             }
166             else
167             {
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.
171                  */
172                 at_ind[f_min] = INT_MAX;
173             }
174         }
175
176         /* Store the bonded end boundaries (at index t) for thread t-1 */
177         for (f = 0; f < ntype; f++)
178         {
179             bt->il_thread_division[ild[f].ftype*(bt->nthreads + 1) + t] = ind[f];
180         }
181     }
182
183     for (f = 0; f < ntype; f++)
184     {
185         assert(ind[f] == ild[f].il->nr);
186     }
187 }
188
189 //! Divides bonded interactions over threads
190 static void divide_bondeds_over_threads(bonded_threading_t *bt,
191                                         const t_idef       *idef)
192 {
193     ilist_data_t ild[F_NRE];
194     int          ntype;
195     int          f;
196
197     assert(bt->nthreads > 0);
198
199     if (F_NRE*(bt->nthreads + 1) > bt->il_thread_division_nalloc)
200     {
201         bt->il_thread_division_nalloc = F_NRE*(bt->nthreads + 1);
202         srenew(bt->il_thread_division, bt->il_thread_division_nalloc);
203     }
204
205     bt->haveBondeds = false;
206     ntype           = 0;
207     for (f = 0; f < F_NRE; f++)
208     {
209         if (!ftype_is_bonded_potential(f))
210         {
211             continue;
212         }
213
214         if (idef->il[f].nr > 0)
215         {
216             bt->haveBondeds = true;
217         }
218
219         if (idef->il[f].nr == 0)
220         {
221             /* No interactions, avoid all the integer math below */
222             int t;
223             for (t = 0; t <= bt->nthreads; t++)
224             {
225                 bt->il_thread_division[f*(bt->nthreads + 1) + t] = 0;
226             }
227         }
228         else if (bt->nthreads <= bt->max_nthread_uniform || f == F_DISRES)
229         {
230             /* On up to 4 threads, load balancing the bonded work
231              * is more important than minimizing the reduction cost.
232              */
233             int nat1, t;
234
235             /* nat1 = 1 + #atoms(ftype) which is the stride use for iatoms */
236             nat1 = 1 + NRAL(f);
237
238             for (t = 0; t <= bt->nthreads; t++)
239             {
240                 int nr_t;
241
242                 /* Divide equally over the threads */
243                 nr_t = (((idef->il[f].nr/nat1)*t)/bt->nthreads)*nat1;
244
245                 if (f == F_DISRES)
246                 {
247                     /* Ensure that distance restraint pairs with the same label
248                      * end up on the same thread.
249                      */
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)
253                     {
254                         nr_t += nat1;
255                     }
256                 }
257
258                 bt->il_thread_division[f*(bt->nthreads + 1) + t] = nr_t;
259             }
260         }
261         else
262         {
263             /* Add this ftype to the list to be distributed */
264             int nat;
265
266             nat              = NRAL(f);
267             ild[ntype].ftype = f;
268             ild[ntype].il    = &idef->il[f];
269             ild[ntype].nat   = nat;
270
271             /* The first index for the thread division is always 0 */
272             bt->il_thread_division[f*(bt->nthreads + 1)] = 0;
273
274             ntype++;
275         }
276     }
277
278     if (ntype > 0)
279     {
280         divide_bondeds_by_locality(bt, ntype, ild);
281     }
282
283     if (debug)
284     {
285         int f;
286
287         fprintf(debug, "Division of bondeds over threads:\n");
288         for (f = 0; f < F_NRE; f++)
289         {
290             if (ftype_is_bonded_potential(f) && idef->il[f].nr > 0)
291             {
292                 int t;
293
294                 fprintf(debug, "%16s", interaction_function[f].name);
295                 for (t = 0; t < bt->nthreads; t++)
296                 {
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])/
300                             (1 + NRAL(f)));
301                 }
302                 fprintf(debug, "\n");
303             }
304         }
305     }
306 }
307
308 //! Construct a reduction mask for which parts (blocks) of the force array are touched on which thread task
309 static void
310 calc_bonded_reduction_mask(int                       natoms,
311                            f_thread_t               *f_thread,
312                            const t_idef             *idef,
313                            int                       thread,
314                            const bonded_threading_t &bondedThreading)
315 {
316     static_assert(BITMASK_SIZE == GMX_OPENMP_MAX_THREADS, "For the error message below we assume these two are equal.");
317
318     if (bondedThreading.nthreads > BITMASK_SIZE)
319     {
320 #pragma omp master
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);
323 #pragma omp barrier
324     }
325     GMX_ASSERT(bondedThreading.nthreads <= BITMASK_SIZE, "We need at least nthreads bits in the mask");
326
327     int nblock = (natoms + reduction_block_size - 1) >> reduction_block_bits;
328
329     if (nblock > f_thread->block_nalloc)
330     {
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);
336     }
337
338     gmx_bitmask_t *mask = f_thread->mask;
339
340     for (int b = 0; b < nblock; b++)
341     {
342         bitmask_clear(&mask[b]);
343     }
344
345     for (int ftype = 0; ftype < F_NRE; ftype++)
346     {
347         if (ftype_is_bonded_potential(ftype))
348         {
349             int nb = idef->il[ftype].nr;
350             if (nb > 0)
351             {
352                 int nat1 = interaction_function[ftype].nratoms + 1;
353
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];
356
357                 for (int i = nb0; i < nb1; i += nat1)
358                 {
359                     for (int a = 1; a < nat1; a++)
360                     {
361                         bitmask_set_bit(&mask[idef->il[ftype].iatoms[i+a] >> reduction_block_bits], thread);
362                     }
363                 }
364             }
365         }
366     }
367
368     /* Make an index of the blocks our thread touches, so we can do fast
369      * force buffer clearing.
370      */
371     f_thread->nblock_used = 0;
372     for (int b = 0; b < nblock; b++)
373     {
374         if (bitmask_is_set(mask[b], thread))
375         {
376             f_thread->block_index[f_thread->nblock_used++] = b;
377         }
378     }
379 }
380
381 void setup_bonded_threading(bonded_threading_t *bt,
382                             int                 numAtoms,
383                             const t_idef       *idef)
384 {
385     int                 ctot = 0;
386
387     assert(bt->nthreads >= 1);
388
389     /* Divide the bonded interaction over the threads */
390     divide_bondeds_over_threads(bt, idef);
391
392     if (!bt->haveBondeds)
393     {
394         /* We don't have bondeds, so there is nothing to reduce */
395         return;
396     }
397
398     /* Determine to which blocks each thread's bonded force calculation
399      * contributes. Store this as a mask for each thread.
400      */
401 #pragma omp parallel for num_threads(bt->nthreads) schedule(static)
402     for (int t = 0; t < bt->nthreads; t++)
403     {
404         try
405         {
406             calc_bonded_reduction_mask(numAtoms, &bt->f_t[t],
407                                        idef, t, *bt);
408         }
409         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
410     }
411
412     /* Reduce the masks over the threads and determine which blocks
413      * we need to reduce over.
414      */
415     int nblock_tot = (numAtoms + reduction_block_size - 1) >> reduction_block_bits;
416     if (nblock_tot > bt->block_nalloc)
417     {
418         bt->block_nalloc = over_alloc_large(nblock_tot);
419         srenew(bt->block_index, bt->block_nalloc);
420         srenew(bt->mask,        bt->block_nalloc);
421     }
422     bt->nblock_used = 0;
423     for (int b = 0; b < nblock_tot; b++)
424     {
425         gmx_bitmask_t *mask = &bt->mask[b];
426
427         /* Generate the union over the threads of the bitmask */
428         bitmask_clear(mask);
429         for (int t = 0; t < bt->nthreads; t++)
430         {
431             bitmask_union(mask, bt->f_t[t].mask[b]);
432         }
433         if (!bitmask_is_zero(*mask))
434         {
435             bt->block_index[bt->nblock_used++] = b;
436         }
437
438         if (debug)
439         {
440             int c = 0;
441             for (int t = 0; t < bt->nthreads; t++)
442             {
443                 if (bitmask_is_set(*mask, t))
444                 {
445                     c++;
446                 }
447             }
448             ctot += c;
449
450             if (gmx_debug_at)
451             {
452 #if BITMASK_SIZE <= 64 //move into bitmask when it is C++
453                 std::string flags = gmx::formatString("%x", *mask);
454 #else
455                 std::string flags = gmx::formatAndJoin(*mask,
456                                                        *mask+BITMASK_ALEN,
457                                                        "", gmx::StringFormatter("%x"));
458 #endif
459                 fprintf(debug, "block %d flags %s count %d\n",
460                         b, flags.c_str(), c);
461             }
462         }
463     }
464     if (debug)
465     {
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);
471     }
472 }
473
474 void tear_down_bonded_threading(bonded_threading_t *bt)
475 {
476     for (int th = 0; th < bt->nthreads; th++)
477     {
478         sfree(bt->f_t[th].fshift);
479         for (int i = 0; i < egNR; i++)
480         {
481             sfree(bt->f_t[th].grpp.ener[i]);
482         }
483     }
484     sfree(bt->f_t);
485     sfree(bt->il_thread_division);
486     sfree(bt);
487 }
488
489 void init_bonded_threading(FILE *fplog, int nenergrp,
490                            struct bonded_threading_t **bt_ptr)
491 {
492     bonded_threading_t *bt;
493
494     snew(bt, 1);
495
496     /* These thread local data structures are used for bondeds only.
497      *
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.
504      */
505     bt->nthreads = gmx_omp_nthreads_get(emntBonded);
506
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++)
510     {
511         try
512         {
513             /* Note that thread 0 uses the global fshift and energy arrays,
514              * but to keep the code simple, we initialize all data here.
515              */
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++)
521             {
522                 snew(bt->f_t[t].grpp.ener[i], bt->f_t[t].grpp.nener);
523             }
524         }
525         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
526     }
527
528     bt->nblock_used  = 0;
529     bt->block_index  = nullptr;
530     bt->mask         = nullptr;
531     bt->block_nalloc = 0;
532
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
535      * and hardware.
536      */
537     const int max_nthread_uniform = 4;
538     char *    ptr;
539
540     if ((ptr = getenv("GMX_BONDED_NTHREAD_UNIFORM")) != nullptr)
541     {
542         sscanf(ptr, "%d", &bt->max_nthread_uniform);
543         if (fplog != nullptr)
544         {
545             fprintf(fplog, "\nMax threads for uniform bonded distribution set to %d by env.var.\n",
546                     bt->max_nthread_uniform);
547         }
548     }
549     else
550     {
551         bt->max_nthread_uniform = max_nthread_uniform;
552     }
553
554     *bt_ptr = bt;
555 }