Avoid using function calls in OpenMP directives
authorErik Lindahl <erik@kth.se>
Tue, 12 Aug 2014 12:15:30 +0000 (14:15 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 30 Sep 2014 11:09:27 +0000 (13:09 +0200)
The direct calls to gmx_omp_nthreads_get() that were
included in some OpenMP pragmas caused memory
corruption and later segfaults on PGI compilers. This
is likely a compiler bug, but we can work around it
by assigning the function return value to a variable
that we use in the pragma.

Such variables are unused when OpenMP is not in use, which might
offend some compiler some time, so adding a gmx_unused attribute is
useful. However, uncrustify needs to be taught about our custom
attributes, which is also done here.

Change-Id: I3b482bdc2401b40a043975ffd4a741f65efd0cfc

admin/uncrustify.cfg
src/gromacs/mdlib/force.c
src/gromacs/mdlib/mdatom.c
src/gromacs/mdlib/minimize.c
src/gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_file_generator/nbnxn_kernel_simd_template.c.pre
src/gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref.c
src/gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.c
src/gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.c
src/gromacs/mdlib/nbnxn_search.c
src/gromacs/mdlib/update.c

index d89351a7e60933183c0ceed7485b270a0e11020d..f7308115b2443d4354ac69a1f903592020f16034 100644 (file)
@@ -1575,3 +1575,12 @@ pp_define_at_level                       = false    # false/true
 # all tokens are separated by any mix of ',' commas, '=' equal signs
 # and whitespace (space, tab)
 #
+
+# Teach uncrustify about the GROMACS attribute aliases that we use
+# to hide compiler differences. This means that declarations like
+#
+# int i, j;
+# int nthreads gmx_unused;
+#
+# does not align i with gmx_unused.
+set ATTRIBUTE gmx_unused gmx_inline gmx_restrict
index 632c2f3a4664ffb00ae72262bcdcf0c6e0dd905a..d4e6445491e36fe2f10bff1cf360d7b6d8ae8050 100644 (file)
@@ -117,9 +117,11 @@ static void reduce_thread_forces(int n, rvec *f,
                                  int nthreads, f_thread_t *f_t)
 {
     int t, i;
+    int nthreads_loop gmx_unused;
 
     /* This reduction can run over any number of threads */
-#pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
+    nthreads_loop = gmx_omp_nthreads_get(emntBonded);
+#pragma omp parallel for num_threads(nthreads_loop) private(t) schedule(static)
     for (i = 0; i < n; i++)
     {
         for (t = 1; t < nthreads; t++)
index 9c0c5943d54f139e40b6629a308b2d59955a2734..b8d51a76da80566855c2088a159bc843d3a78e4e 100644 (file)
@@ -121,6 +121,7 @@ void atoms2md(gmx_mtop_t *mtop, t_inputrec *ir,
     t_grpopts            *opts;
     gmx_groups_t         *groups;
     gmx_molblock_t       *molblock;
+    int                   nthreads gmx_unused;
 
     bLJPME = EVDW_PME(ir->vdwtype);
 
@@ -233,7 +234,8 @@ void atoms2md(gmx_mtop_t *mtop, t_inputrec *ir,
 
     alook = gmx_mtop_atomlookup_init(mtop);
 
-#pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntDefault)) schedule(static)
+    nthreads = gmx_omp_nthreads_get(emntDefault);
+#pragma omp parallel for num_threads(nthreads) schedule(static)
     for (i = 0; i < md->nr; i++)
     {
         int      g, ag, molb;
index 3cc1bbc8a08c6a9d65f6d324369dc49cc6addcd3..69008f53fabcb3388951c01d7554d479c67ccdca 100644 (file)
@@ -550,6 +550,7 @@ static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
     int      start, end;
     rvec    *x1, *x2;
     real     dvdl_constr;
+    int      nthreads gmx_unused;
 
     s1 = &ems1->s;
     s2 = &ems2->s;
@@ -587,7 +588,8 @@ static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
     x1 = s1->x;
     x2 = s2->x;
 
-#pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
+    nthreads = gmx_omp_nthreads_get(emntUpdate);
+#pragma omp parallel num_threads(nthreads)
     {
         int gf, i, m;
 
index a3ede7190f921840e851d48b81fb5554c8bb5365..2aa4fa9fb59bf5db2dab19d027552d4707556867 100644 (file)
@@ -138,6 +138,7 @@ void
     nbnxn_pairlist_t **nbl;
     int                coulkt, vdwkt = 0;
     int                nb;
+    int                nthreads gmx_unused;
 
     nnbl = nbl_list->nnbl;
     nbl  = nbl_list->nbl;
@@ -209,7 +210,8 @@ void
         gmx_incons("Unsupported VdW interaction type");
     }}
 
-#pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
+    nthreads = gmx_omp_nthreads_get(emntNonbonded);
+#pragma omp parallel for schedule(static) num_threads(nthreads)
     for (nb = 0; nb < nnbl; nb++)
     {{
         nbnxn_atomdata_output_t *out;
index 0be4be644d987e9eecf22b6a64f4457644cfe755..a10465296569d932b4b90a2b34ca131d98c2fe58 100644 (file)
@@ -182,6 +182,7 @@ nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
     int                coult;
     int                vdwt;
     int                nb;
+    int                nthreads gmx_unused;
 
     nnbl = nbl_list->nnbl;
     nbl  = nbl_list->nbl;
@@ -239,7 +240,8 @@ nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
         gmx_incons("Unsupported vdwtype in nbnxn reference kernel");
     }
 
-#pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
+    nthreads = gmx_omp_nthreads_get(emntNonbonded);
+#pragma omp parallel for schedule(static) num_threads(nthreads)
     for (nb = 0; nb < nnbl; nb++)
     {
         nbnxn_atomdata_output_t *out;
index 45f09bc4621d64d900ba6c9fc66effead4372743..35181d43eb17c584fa2af44b62f20aef15443875 100644 (file)
@@ -273,6 +273,7 @@ nbnxn_kernel_simd_2xnn(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
     nbnxn_pairlist_t **nbl;
     int                coulkt, vdwkt = 0;
     int                nb;
+    int                nthreads gmx_unused;
 
     nnbl = nbl_list->nnbl;
     nbl  = nbl_list->nbl;
@@ -344,7 +345,8 @@ nbnxn_kernel_simd_2xnn(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
         gmx_incons("Unsupported VdW interaction type");
     }
 
-#pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
+    nthreads = gmx_omp_nthreads_get(emntNonbonded);
+#pragma omp parallel for schedule(static) num_threads(nthreads)
     for (nb = 0; nb < nnbl; nb++)
     {
         nbnxn_atomdata_output_t *out;
index 0684149bad85f093838e42c941c8593d7aa05567..5974873a070147da2fd1723c8909d00449eb3dfe 100644 (file)
@@ -272,6 +272,7 @@ nbnxn_kernel_simd_4xn(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
     nbnxn_pairlist_t **nbl;
     int                coulkt, vdwkt = 0;
     int                nb;
+    int                nthreads gmx_unused;
 
     nnbl = nbl_list->nnbl;
     nbl  = nbl_list->nbl;
@@ -343,7 +344,8 @@ nbnxn_kernel_simd_4xn(nbnxn_pairlist_set_t      gmx_unused *nbl_list,
         gmx_incons("Unsupported VdW interaction type");
     }
 
-#pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
+    nthreads = gmx_omp_nthreads_get(emntNonbonded);
+#pragma omp parallel for schedule(static) num_threads(nthreads)
     for (nb = 0; nb < nnbl; nb++)
     {
         nbnxn_atomdata_output_t *out;
index dc0b02013f645ed2f44acb66c6509f23e5350cb5..3114bebc22b048559369dede69c9774cf4d8f21a 100644 (file)
@@ -1931,6 +1931,7 @@ void nbnxn_grid_add_simple(nbnxn_search_t    nbs,
     float        *bbcz;
     nbnxn_bb_t   *bb;
     int           ncd, sc;
+    int           nthreads gmx_unused;
 
     grid = &nbs->grid[0];
 
@@ -1957,7 +1958,8 @@ void nbnxn_grid_add_simple(nbnxn_search_t    nbs,
     bbcz = grid->bbcz_simple;
     bb   = grid->bb_simple;
 
-#pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
+    nthreads = gmx_omp_nthreads_get(emntPairsearch);
+#pragma omp parallel for num_threads(nthreads) schedule(static)
     for (sc = 0; sc < grid->nc; sc++)
     {
         int c, tx, na;
@@ -4473,6 +4475,7 @@ static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
 {
     int nsci, ncj4, nexcl;
     int n, i;
+    int nthreads gmx_unused;
 
     if (nblc->bSimple)
     {
@@ -4513,7 +4516,8 @@ static void combine_nblists(int nnbl, nbnxn_pairlist_t **nbl,
     /* Each thread should copy its own data to the combined arrays,
      * as otherwise data will go back and forth between different caches.
      */
-#pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntPairsearch)) schedule(static)
+    nthreads = gmx_omp_nthreads_get(emntPairsearch);
+#pragma omp parallel for num_threads(nthreads) schedule(static)
     for (n = 0; n < nnbl; n++)
     {
         int                     sci_offset;
index 787135f7158cd12b23bac96a7446d7de8d8e44a7..93b7f59c1ec4dc7faffce100b5b1f1f51920c478 100644 (file)
@@ -1771,7 +1771,9 @@ void update_constraints(FILE             *fplog,
         }
         else
         {
-#pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntUpdate)) schedule(static)
+            nth = gmx_omp_nthreads_get(emntUpdate);
+
+#pragma omp parallel for num_threads(nth) schedule(static)
             for (i = start; i < nrend; i++)
             {
                 copy_rvec(upd->xp[i], state->x[i]);