Add InteractionDefinitions
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /*! \internal \file
39  * \brief This file defines functions for managing threading of listed
40  * interactions.
41  *
42  * \author Mark Abraham <mark.j.abraham@gmail.com>
43  * \ingroup module_listed_forces
44  */
45 #include "gmxpre.h"
46
47 #include "manage_threading.h"
48
49 #include "config.h"
50
51 #include <cassert>
52 #include <cinttypes>
53 #include <climits>
54 #include <cstdlib>
55
56 #include <algorithm>
57 #include <string>
58
59 #include "gromacs/listed_forces/gpubonded.h"
60 #include "gromacs/mdlib/gmx_omp_nthreads.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/topology/ifunc.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/smalloc.h"
67
68 #include "listed_internal.h"
69 #include "utilities.h"
70
71 /*! \brief struct for passing all data required for a function type */
72 typedef struct
73 {
74     const InteractionList* il;    /**< pointer to t_ilist entry corresponding to ftype */
75     int                    ftype; /**< the function type index */
76     int                    nat;   /**< nr of atoms involved in a single ftype interaction */
77 } ilist_data_t;
78
79 /*! \brief Divides listed interactions over threads
80  *
81  * This routine attempts to divide all interactions of the numType bondeds
82  * types stored in ild over the threads such that each thread has roughly
83  * equal load and different threads avoid touching the same atoms as much
84  * as possible.
85  */
86 static void divide_bondeds_by_locality(bonded_threading_t* bt, int numType, const ilist_data_t* ild)
87 {
88     int nat_tot, nat_sum;
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 */
91     int f, t;
92
93     assert(numType <= F_NRE);
94
95     nat_tot = 0;
96     for (f = 0; f < numType; f++)
97     {
98         /* Sum #bondeds*#atoms_per_bond over all bonded types */
99         nat_tot += ild[f].il->size() / (ild[f].nat + 1) * ild[f].nat;
100         /* The start bound for thread 0 is 0 for all interactions */
101         ind[f] = 0;
102         /* Initialize the next atom index array */
103         assert(!ild[f].il->empty());
104         at_ind[f] = ild[f].il->iatoms[1];
105     }
106
107     nat_sum = 0;
108     /* Loop over the end bounds of the nthreads threads to determine
109      * which interactions threads 0 to nthreads shall calculate.
110      *
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.
118      */
119     for (t = 1; t <= bt->nthreads; t++)
120     {
121         int nat_thread;
122
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.
129          */
130         nat_thread = (nat_tot * t) / bt->nthreads;
131
132         while (nat_sum < nat_thread)
133         {
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
143              * be higher.
144              */
145             int f_min;
146
147             /* Find out which of the types has the lowest atom index */
148             f_min = 0;
149             for (f = 1; f < numType; f++)
150             {
151                 if (at_ind[f] < at_ind[f_min])
152                 {
153                     f_min = f;
154                 }
155             }
156             assert(f_min >= 0 && f_min < numType);
157
158             /* Assign the interaction with the lowest atom index (of type
159              * index f_min) to thread t-1 by increasing ind.
160              */
161             ind[f_min] += ild[f_min].nat + 1;
162             nat_sum += ild[f_min].nat;
163
164             /* Update the first unassigned atom index for this type */
165             if (ind[f_min] < ild[f_min].il->size())
166             {
167                 at_ind[f_min] = ild[f_min].il->iatoms[ind[f_min] + 1];
168             }
169             else
170             {
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.
174                  */
175                 at_ind[f_min] = INT_MAX;
176             }
177         }
178
179         /* Store the bonded end boundaries (at index t) for thread t-1 */
180         for (f = 0; f < numType; f++)
181         {
182             bt->workDivision.setBound(ild[f].ftype, t, ind[f]);
183         }
184     }
185
186     for (f = 0; f < numType; f++)
187     {
188         assert(ind[f] == ild[f].il->size());
189     }
190 }
191
192 //! Return whether function type \p ftype in \p idef has perturbed interactions
193 static bool ftypeHasPerturbedEntries(const InteractionDefinitions& idef, int ftype)
194 {
195     GMX_ASSERT(idef.ilsort == ilsortNO_FE || idef.ilsort == ilsortFE_SORTED,
196                "Perturbed interations should be sorted here");
197
198     const InteractionList& ilist = idef.il[ftype];
199
200     return (idef.ilsort != ilsortNO_FE && idef.numNonperturbedInteractions[ftype] != ilist.size());
201 }
202
203 //! Divides bonded interactions over threads and GPU
204 static void divide_bondeds_over_threads(bonded_threading_t*           bt,
205                                         bool                          useGpuForBondeds,
206                                         const InteractionDefinitions& idef)
207 {
208     ilist_data_t ild[F_NRE];
209
210     GMX_ASSERT(bt->nthreads > 0, "Must have positive number of threads");
211     const int numThreads = bt->nthreads;
212
213     gmx::ArrayRef<const t_iparams> iparams = idef.iparams;
214
215     bt->haveBondeds      = false;
216     int    numType       = 0;
217     size_t fTypeGpuIndex = 0;
218     for (int fType = 0; fType < F_NRE; fType++)
219     {
220         if (!ftype_is_bonded_potential(fType))
221         {
222             continue;
223         }
224
225         const InteractionList& il                     = idef.il[fType];
226         int                    nrToAssignToCpuThreads = il.size();
227
228         if (useGpuForBondeds && fTypeGpuIndex < gmx::fTypesOnGpu.size()
229             && gmx::fTypesOnGpu[fTypeGpuIndex] == fType)
230         {
231             fTypeGpuIndex++;
232
233             /* Perturbation is not implemented in the GPU bonded kernels.
234              * But instead of doing all on the CPU, we could do only
235              * the actually perturbed interactions on the CPU.
236              */
237             if (!ftypeHasPerturbedEntries(idef, fType))
238             {
239                 /* We will assign this interaction type to the GPU */
240                 nrToAssignToCpuThreads = 0;
241             }
242         }
243
244         if (nrToAssignToCpuThreads > 0)
245         {
246             bt->haveBondeds = true;
247         }
248
249         if (nrToAssignToCpuThreads == 0)
250         {
251             /* No interactions, avoid all the integer math below */
252             for (int t = 0; t <= numThreads; t++)
253             {
254                 bt->workDivision.setBound(fType, t, 0);
255             }
256         }
257         else if (numThreads <= bt->max_nthread_uniform || fType == F_DISRES)
258         {
259             /* On up to 4 threads, load balancing the bonded work
260              * is more important than minimizing the reduction cost.
261              */
262
263             const int stride = 1 + NRAL(fType);
264
265             for (int t = 0; t <= numThreads; t++)
266             {
267                 /* Divide equally over the threads */
268                 int nr_t = (((nrToAssignToCpuThreads / stride) * t) / numThreads) * stride;
269
270                 if (fType == F_DISRES)
271                 {
272                     /* Ensure that distance restraint pairs with the same label
273                      * end up on the same thread.
274                      */
275                     while (nr_t > 0 && nr_t < nrToAssignToCpuThreads
276                            && iparams[il.iatoms[nr_t]].disres.label
277                                       == iparams[il.iatoms[nr_t - stride]].disres.label)
278                     {
279                         nr_t += stride;
280                     }
281                 }
282
283                 bt->workDivision.setBound(fType, t, nr_t);
284             }
285         }
286         else
287         {
288             /* Add this fType to the list to be distributed */
289             int nat            = NRAL(fType);
290             ild[numType].ftype = fType;
291             ild[numType].il    = &il;
292             ild[numType].nat   = nat;
293
294             /* The first index for the thread division is always 0 */
295             bt->workDivision.setBound(fType, 0, 0);
296
297             numType++;
298         }
299     }
300
301     if (numType > 0)
302     {
303         divide_bondeds_by_locality(bt, numType, ild);
304     }
305
306     if (debug)
307     {
308         int f;
309
310         fprintf(debug, "Division of bondeds over threads:\n");
311         for (f = 0; f < F_NRE; f++)
312         {
313             if (ftype_is_bonded_potential(f) && !idef.il[f].empty())
314             {
315                 int t;
316
317                 fprintf(debug, "%16s", interaction_function[f].name);
318                 for (t = 0; t < numThreads; t++)
319                 {
320                     fprintf(debug, " %4d",
321                             (bt->workDivision.bound(f, t + 1) - bt->workDivision.bound(f, t))
322                                     / (1 + NRAL(f)));
323                 }
324                 fprintf(debug, "\n");
325             }
326         }
327     }
328 }
329
330 //! Construct a reduction mask for which parts (blocks) of the force array are touched on which thread task
331 static void calc_bonded_reduction_mask(int                           natoms,
332                                        f_thread_t*                   f_thread,
333                                        const InteractionDefinitions& idef,
334                                        int                           thread,
335                                        const bonded_threading_t&     bondedThreading)
336 {
337     static_assert(BITMASK_SIZE == GMX_OPENMP_MAX_THREADS,
338                   "For the error message below we assume these two are equal.");
339
340     if (bondedThreading.nthreads > BITMASK_SIZE)
341     {
342 #pragma omp master
343         gmx_fatal(FARGS,
344                   "You are using %d OpenMP threads, which is larger than GMX_OPENMP_MAX_THREADS "
345                   "(%d). Decrease the number of OpenMP threads or rebuild GROMACS with a larger "
346                   "value for GMX_OPENMP_MAX_THREADS passed to CMake.",
347                   bondedThreading.nthreads, GMX_OPENMP_MAX_THREADS);
348 #pragma omp barrier
349     }
350     GMX_ASSERT(bondedThreading.nthreads <= BITMASK_SIZE,
351                "We need at least nthreads bits in the mask");
352
353     int nblock = (natoms + reduction_block_size - 1) >> reduction_block_bits;
354
355     if (nblock > f_thread->block_nalloc)
356     {
357         f_thread->block_nalloc = over_alloc_large(nblock);
358         srenew(f_thread->mask, f_thread->block_nalloc);
359         srenew(f_thread->block_index, f_thread->block_nalloc);
360         // NOTE: It seems f_thread->f does not need to be aligned
361         sfree_aligned(f_thread->f);
362         snew_aligned(f_thread->f, f_thread->block_nalloc * reduction_block_size, 128);
363     }
364
365     gmx_bitmask_t* mask = f_thread->mask;
366
367     for (int b = 0; b < nblock; b++)
368     {
369         bitmask_clear(&mask[b]);
370     }
371
372     for (int ftype = 0; ftype < F_NRE; ftype++)
373     {
374         if (ftype_is_bonded_potential(ftype))
375         {
376             int nb = idef.il[ftype].size();
377             if (nb > 0)
378             {
379                 int nat1 = interaction_function[ftype].nratoms + 1;
380
381                 int nb0 = bondedThreading.workDivision.bound(ftype, thread);
382                 int nb1 = bondedThreading.workDivision.bound(ftype, thread + 1);
383
384                 for (int i = nb0; i < nb1; i += nat1)
385                 {
386                     for (int a = 1; a < nat1; a++)
387                     {
388                         bitmask_set_bit(&mask[idef.il[ftype].iatoms[i + a] >> reduction_block_bits], thread);
389                     }
390                 }
391             }
392         }
393     }
394
395     /* Make an index of the blocks our thread touches, so we can do fast
396      * force buffer clearing.
397      */
398     f_thread->nblock_used = 0;
399     for (int b = 0; b < nblock; b++)
400     {
401         if (bitmask_is_set(mask[b], thread))
402         {
403             f_thread->block_index[f_thread->nblock_used++] = b;
404         }
405     }
406 }
407
408 void setup_bonded_threading(bonded_threading_t*           bt,
409                             int                           numAtoms,
410                             bool                          useGpuForBondeds,
411                             const InteractionDefinitions& idef)
412 {
413     int ctot = 0;
414
415     assert(bt->nthreads >= 1);
416
417     /* Divide the bonded interaction over the threads */
418     divide_bondeds_over_threads(bt, useGpuForBondeds, idef);
419
420     if (!bt->haveBondeds)
421     {
422         /* We don't have bondeds, so there is nothing to reduce */
423         return;
424     }
425
426     /* Determine to which blocks each thread's bonded force calculation
427      * contributes. Store this as a mask for each thread.
428      */
429 #pragma omp parallel for num_threads(bt->nthreads) schedule(static)
430     for (int t = 0; t < bt->nthreads; t++)
431     {
432         try
433         {
434             calc_bonded_reduction_mask(numAtoms, bt->f_t[t].get(), idef, t, *bt);
435         }
436         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
437     }
438
439     /* Reduce the masks over the threads and determine which blocks
440      * we need to reduce over.
441      */
442     int nblock_tot = (numAtoms + reduction_block_size - 1) >> reduction_block_bits;
443     /* Ensure we have sufficient space for all blocks */
444     if (static_cast<size_t>(nblock_tot) > bt->block_index.size())
445     {
446         bt->block_index.resize(nblock_tot);
447     }
448     if (static_cast<size_t>(nblock_tot) > bt->mask.size())
449     {
450         bt->mask.resize(nblock_tot);
451     }
452     bt->nblock_used = 0;
453     for (int b = 0; b < nblock_tot; b++)
454     {
455         gmx_bitmask_t* mask = &bt->mask[b];
456
457         /* Generate the union over the threads of the bitmask */
458         bitmask_clear(mask);
459         for (int t = 0; t < bt->nthreads; t++)
460         {
461             bitmask_union(mask, bt->f_t[t]->mask[b]);
462         }
463         if (!bitmask_is_zero(*mask))
464         {
465             bt->block_index[bt->nblock_used++] = b;
466         }
467
468         if (debug)
469         {
470             int c = 0;
471             for (int t = 0; t < bt->nthreads; t++)
472             {
473                 if (bitmask_is_set(*mask, t))
474                 {
475                     c++;
476                 }
477             }
478             ctot += c;
479
480             if (gmx_debug_at)
481             {
482                 fprintf(debug, "block %d flags %s count %d\n", b, to_hex_string(*mask).c_str(), c);
483             }
484         }
485     }
486     if (debug)
487     {
488         fprintf(debug, "Number of %d atom blocks to reduce: %d\n", reduction_block_size, bt->nblock_used);
489         fprintf(debug, "Reduction density %.2f for touched blocks only %.2f\n",
490                 ctot * reduction_block_size / static_cast<double>(numAtoms),
491                 ctot / static_cast<double>(bt->nblock_used));
492     }
493 }
494
495 void tear_down_bonded_threading(bonded_threading_t* bt)
496 {
497     delete bt;
498 }
499
500 f_thread_t::f_thread_t(int numEnergyGroups) : grpp(numEnergyGroups)
501 {
502     snew(fshift, SHIFTS);
503 }
504
505 f_thread_t::~f_thread_t()
506 {
507     sfree(mask);
508     sfree(fshift);
509     sfree(block_index);
510     sfree_aligned(f);
511 }
512
513 bonded_threading_t::bonded_threading_t(const int numThreads, const int numEnergyGroups) :
514     nthreads(numThreads),
515     nblock_used(0),
516     haveBondeds(false),
517     workDivision(nthreads),
518     foreignLambdaWorkDivision(1)
519 {
520     f_t.resize(numThreads);
521 #pragma omp parallel for num_threads(nthreads) schedule(static)
522     for (int t = 0; t < nthreads; t++)
523     {
524         try
525         {
526             /* Note that thread 0 uses the global fshift and energy arrays,
527              * but to keep the code simple, we initialize all data here.
528              */
529             f_t[t] = std::make_unique<f_thread_t>(numEnergyGroups);
530         }
531         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
532     }
533 }
534
535 bonded_threading_t* init_bonded_threading(FILE* fplog, const int nenergrp)
536 {
537     /* These thread local data structures are used for bondeds only.
538      *
539      * Note that we also use there structures when running single-threaded.
540      * This is because the bonded force buffer uses type rvec4, whereas
541      * the normal force buffer is uses type rvec. This leads to a little
542      * reduction overhead, but the speed gain in the bonded calculations
543      * of doing transposeScatterIncr/DecrU with aligment 4 instead of 3
544      * is much larger than the reduction overhead.
545      */
546     bonded_threading_t* bt = new bonded_threading_t(gmx_omp_nthreads_get(emntBonded), nenergrp);
547
548     /* The optimal value after which to switch from uniform to localized
549      * bonded interaction distribution is 3, 4 or 5 depending on the system
550      * and hardware.
551      */
552     const int max_nthread_uniform = 4;
553     char*     ptr;
554
555     if ((ptr = getenv("GMX_BONDED_NTHREAD_UNIFORM")) != nullptr)
556     {
557         sscanf(ptr, "%d", &bt->max_nthread_uniform);
558         if (fplog != nullptr)
559         {
560             fprintf(fplog, "\nMax threads for uniform bonded distribution set to %d by env.var.\n",
561                     bt->max_nthread_uniform);
562         }
563     }
564     else
565     {
566         bt->max_nthread_uniform = max_nthread_uniform;
567     }
568
569     return bt;
570 }