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, 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.
39 #include "gromacs/legacyheaders/bonded-threading.h"
45 #include "gromacs/bonded/bonded.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
49 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
56 idef->nthreads = nthreads;
58 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
60 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
61 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
64 for (ftype = 0; ftype < F_NRE; ftype++)
66 if (ftype_is_bonded_potential(ftype))
68 nat1 = interaction_function[ftype].nratoms + 1;
70 for (t = 0; t <= nthreads; t++)
72 /* Divide the interactions equally over the threads.
73 * When the different types of bonded interactions
74 * are distributed roughly equally over the threads,
75 * this should lead to well localized output into
76 * the force buffer on each thread.
77 * If this is not the case, a more advanced scheme
78 * (not implemented yet) will do better.
80 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
82 /* Ensure that distance restraint pairs with the same label
83 * end up on the same thread.
84 * This is slighlty tricky code, since the next for iteration
85 * may have an initial il_nr_thread lower than the final value
86 * in the previous iteration, but this will anyhow be increased
87 * to the approriate value again by this while loop.
89 while (ftype == F_DISRES &&
91 il_nr_thread < idef->il[ftype].nr &&
92 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
93 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
98 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
105 calc_bonded_reduction_mask(const t_idef *idef,
110 int ftype, nb, nat1, nb0, nb1, i, a;
114 for (ftype = 0; ftype < F_NRE; ftype++)
116 if (ftype_is_bonded_potential(ftype))
118 nb = idef->il[ftype].nr;
121 nat1 = interaction_function[ftype].nratoms + 1;
123 /* Divide this interaction equally over the threads.
124 * This is not stored: should match division in calc_bonds.
126 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
127 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
129 for (i = nb0; i < nb1; i += nat1)
131 for (a = 1; a < nat1; a++)
133 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
143 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
145 #define MAX_BLOCK_BITS 32
149 assert(fr->nthreads >= 1);
151 /* Divide the bonded interaction over the threads */
152 divide_bondeds_over_threads(idef, fr->nthreads);
154 if (fr->nthreads == 1)
161 /* We divide the force array in a maximum of 32 blocks.
162 * Minimum force block reduction size is 2^6=64.
165 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
171 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
175 /* Determine to which blocks each thread's bonded force calculation
176 * contributes. Store this is a mask for each thread.
178 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
179 for (t = 1; t < fr->nthreads; t++)
181 fr->f_t[t].red_mask =
182 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
185 /* Determine the maximum number of blocks we need to reduce over */
188 for (t = 0; t < fr->nthreads; t++)
191 for (b = 0; b < MAX_BLOCK_BITS; b++)
193 if (fr->f_t[t].red_mask & (1U<<b))
195 fr->red_nblock = std::max(fr->red_nblock, b+1);
201 fprintf(debug, "thread %d flags %x count %d\n",
202 t, fr->f_t[t].red_mask, c);
208 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
209 fr->red_nblock, 1<<fr->red_ashift);
210 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
211 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
212 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));