Clean up thread teardown
authorKevin Boyd <kevin.boyd@uconn.edu>
Mon, 7 Jan 2019 22:41:24 +0000 (17:41 -0500)
committerPaul Bauer <paul.bauer.q@gmail.com>
Wed, 9 Jan 2019 05:58:32 +0000 (06:58 +0100)
Struct freeing had been incomplete

Change-Id: If48f495f64bf880714b26241aa5571f5d487521b

src/gromacs/listed-forces/listed-internal.h
src/gromacs/listed-forces/manage-threading.cpp

index 87f7197255306523fb50a1aff5b2eaea0b4493be..01b11aa255b53c35c5c20d3fc5450544a55682b0 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014,2015,2016,2018, by the GROMACS development team, led by
+ * Copyright (c) 2014,2015,2016,2018,2019, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -78,7 +78,7 @@ struct bonded_threading_t
 {
     /* Thread local force and energy data */
     int            nthreads;     /**< Number of threads to be used for bondeds */
-    f_thread_t    *f_t;          /**< Force/enegry data per thread, size nthreads */
+    f_thread_t    *f_t;          /**< Force/energy data per thread, size nthreads */
     int            nblock_used;  /**< The number of force blocks to reduce */
     int           *block_index;  /**< Index of size nblock_used into mask */
     gmx_bitmask_t *mask;         /**< Mask array, one element corresponds to a block of reduction_block_size atoms of the force array, bit corresponding to thread indices set if a thread writes to that block */
index 3a2e54e3c32d399df7ef3b3c7c3d3e47b7859a67..c1fc232fc8f2f4ef6e377a5a9bca3c0ca585b782 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -497,7 +497,10 @@ void tear_down_bonded_threading(bonded_threading_t *bt)
 {
     for (int th = 0; th < bt->nthreads; th++)
     {
+        sfree(bt->f_t[th].mask);
         sfree(bt->f_t[th].fshift);
+        sfree(bt->f_t[th].block_index);
+        sfree_aligned(bt->f_t[th].f);
         for (int i = 0; i < egNR; i++)
         {
             sfree(bt->f_t[th].grpp.ener[i]);