Moving put_atoms_in_box_omp() to pbc.h
authorPrashanth Kanduri <kanduri@cscs.ch>
Tue, 19 Mar 2019 15:37:46 +0000 (16:37 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 27 Mar 2019 13:33:27 +0000 (14:33 +0100)
This is another patch in the cleaning efforts of sim_util.
It finally removes sim_util.h

The idea is to only have functions relevant to the force
schedules there so that it becomes easy to move it into
its own module with minimal merging pains.

Related: #2574

Change-Id: Ib0c2dbc21bd31ee272888d3fa25a3c0ce65b5478

src/gromacs/gmxpreprocess/grompp.cpp
src/gromacs/mdlib/forcerec.cpp
src/gromacs/mdlib/shellfc.cpp
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdlib/sim_util.h [deleted file]
src/gromacs/mdrun/md.cpp
src/gromacs/mdrun/tpi.cpp
src/gromacs/pbcutil/pbc.cpp
src/gromacs/pbcutil/pbc.h
src/gromacs/swap/swapcoords.cpp

index 8b01bd5f25ac3c71769bdd1288ef0db50cac5edb..189743bdc13841cec63753137cc35986fd376798 100644 (file)
@@ -80,7 +80,6 @@
 #include "gromacs/mdlib/constr.h"
 #include "gromacs/mdlib/perf_est.h"
 #include "gromacs/mdlib/qmmm.h"
-#include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/vsite.h"
 #include "gromacs/mdrunutility/mdmodules.h"
 #include "gromacs/mdtypes/inputrec.h"
index e852525ebd7201a6bcd5ef4f3ce4f2a254f90804..1d42ce69b2296f31094af903df90cb51e50e5d8b 100644 (file)
@@ -73,7 +73,6 @@
 #include "gromacs/mdlib/ns.h"
 #include "gromacs/mdlib/qmmm.h"
 #include "gromacs/mdlib/rf_util.h"
-#include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/wall.h"
 #include "gromacs/mdtypes/commrec.h"
 #include "gromacs/mdtypes/fcdata.h"
index 76caff90036a0963955a78c90fe3523fffa12530..9fa062db9c6c269e0e286bc4f17d9b94b5d38ac1 100644 (file)
@@ -58,8 +58,8 @@
 #include "gromacs/mdlib/constr.h"
 #include "gromacs/mdlib/force.h"
 #include "gromacs/mdlib/force_flags.h"
+#include "gromacs/mdlib/gmx_omp_nthreads.h"
 #include "gromacs/mdlib/mdatoms.h"
-#include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/vsite.h"
 #include "gromacs/mdtypes/commrec.h"
 #include "gromacs/mdtypes/enerdata.h"
@@ -1059,7 +1059,7 @@ void relax_shell_flexcon(FILE                                     *fplog,
         if (inputrec->cutoff_scheme == ecutsVERLET)
         {
             auto xRef = state->x.arrayRefWithPadding().paddedArrayRef();
-            put_atoms_in_box_omp(fr->ePBC, state->box, xRef.subArray(0, md->homenr));
+            put_atoms_in_box_omp(fr->ePBC, state->box, xRef.subArray(0, md->homenr), gmx_omp_nthreads_get(emntDefault));
         }
         else
         {
index 1c3fc508d042c61ee7aef75af915bfd622c494c0..843b1db72d068f9bdbafda7aacd43f200aa98b92 100644 (file)
@@ -36,8 +36,6 @@
  */
 #include "gmxpre.h"
 
-#include "sim_util.h"
-
 #include "config.h"
 
 #include <cmath>
@@ -892,7 +890,7 @@ static void do_force_cutsVERLET(FILE *fplog,
 
         if (bCalcCGCM)
         {
-            put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr));
+            put_atoms_in_box_omp(fr->ePBC, box, x.unpaddedArrayRef().subArray(0, homenr), gmx_omp_nthreads_get(emntDefault));
             inc_nrnb(nrnb, eNR_SHIFTX, homenr);
         }
         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
@@ -1869,22 +1867,3 @@ void do_force(FILE                                     *fplog,
      */
     ddBalanceRegionHandler.openBeforeForceComputationCpu(DdAllowBalanceRegionReopen::no);
 }
-
-void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x)
-{
-    int t, nth;
-    nth = gmx_omp_nthreads_get(emntDefault);
-
-#pragma omp parallel for num_threads(nth) schedule(static)
-    for (t = 0; t < nth; t++)
-    {
-        try
-        {
-            size_t natoms = x.size();
-            size_t offset = (natoms*t    )/nth;
-            size_t len    = (natoms*(t + 1))/nth - offset;
-            put_atoms_in_box(ePBC, box, x.subArray(offset, len));
-        }
-        GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
-    }
-}
diff --git a/src/gromacs/mdlib/sim_util.h b/src/gromacs/mdlib/sim_util.h
deleted file mode 100644 (file)
index 036aac4..0000000
+++ /dev/null
@@ -1,54 +0,0 @@
-/*
- * This file is part of the GROMACS molecular simulation package.
- *
- * 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,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.
- *
- * GROMACS is free software; you can redistribute it and/or
- * modify it under the terms of the GNU Lesser General Public License
- * as published by the Free Software Foundation; either version 2.1
- * of the License, or (at your option) any later version.
- *
- * GROMACS is distributed in the hope that it will be useful,
- * but WITHOUT ANY WARRANTY; without even the implied warranty of
- * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
- * Lesser General Public License for more details.
- *
- * You should have received a copy of the GNU Lesser General Public
- * License along with GROMACS; if not, see
- * http://www.gnu.org/licenses, or write to the Free Software Foundation,
- * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
- *
- * If you want to redistribute modifications to GROMACS, please
- * consider that scientific software is very special. Version
- * control is crucial - bugs must be traceable. We will be happy to
- * consider code for inclusion in the official distribution, but
- * derived work must not be called official GROMACS. Details are found
- * in the README & COPYING files - if they are missing, get the
- * official version at http://www.gromacs.org.
- *
- * To help us fund GROMACS development, we humbly ask that you cite
- * the research papers on the package. Check out http://www.gromacs.org.
- */
-#ifndef GMX_MDLIB_SIM_UTIL_H
-#define GMX_MDLIB_SIM_UTIL_H
-
-#include "gromacs/math/vectypes.h"
-#include "gromacs/utility/arrayref.h"
-
-/*! \brief Parallellizes put_atoms_in_box()
- *
- * This wrapper function around put_atoms_in_box() with the ugly manual
- * workload splitting is needed to avoid silently introducing multithreading
- * in tools.
- * \param[in]    ePBC   The pbc type
- * \param[in]    box    The simulation box
- * \param[inout] x      The coordinates of the atoms
- */
-void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x);
-
-#endif
index 338ff42958ad853d1be2ae7b2eb210504f8d889a..5e1db3704fc98b51ea4a078ee82540d02a62d1ef 100644 (file)
@@ -90,7 +90,6 @@
 #include "gromacs/mdlib/resethandler.h"
 #include "gromacs/mdlib/shellfc.h"
 #include "gromacs/mdlib/sighandler.h"
-#include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/simulationsignal.h"
 #include "gromacs/mdlib/stat.h"
 #include "gromacs/mdlib/stophandler.h"
index 227ce96fd7529f925874c7074f6bf5ba48d550c7..cb75faf26a31189a79d144b4607e853ec9218456 100644 (file)
@@ -72,7 +72,6 @@
 #include "gromacs/mdlib/force_flags.h"
 #include "gromacs/mdlib/mdatoms.h"
 #include "gromacs/mdlib/ns.h"
-#include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdlib/tgroup.h"
 #include "gromacs/mdlib/update.h"
 #include "gromacs/mdlib/vsite.h"
index b77448fc5f6c978c1356d83465d26cfbc0fb046a..0397683415186f86a1bc23231bff768680814190 100644 (file)
@@ -1471,6 +1471,24 @@ void put_atoms_in_box(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x)
     }
 }
 
+void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth)
+{
+    int t;
+
+#pragma omp parallel for num_threads(nth) schedule(static)
+    for (t = 0; t < nth; t++)
+    {
+        try
+        {
+            size_t natoms = x.size();
+            size_t offset = (natoms*t    )/nth;
+            size_t len    = (natoms*(t + 1))/nth - offset;
+            put_atoms_in_box(ePBC, box, x.subArray(offset, len));
+        }
+        GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
+    }
+}
+
 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box,
                                      gmx::ArrayRef<gmx::RVec> x)
 {
index ad8d41ff05d8ead686d2a8f25d5c97489cb4d61f..dc70f22994e26a4cfc44d2f1d861dbb794f0d959 100644 (file)
@@ -306,12 +306,24 @@ int *compact_unitcell_edges();
  * These routines puts ONE or ALL atoms in the box, not caring
  * about charge groups!
  * Also works for triclinic cells.
- * \param[in]    ePBC   The pbc type
- * \param[in]    box    The simulation box
+ * \param[in]     ePBC   The pbc type
+ * \param[in]     box    The simulation box
  * \param[in,out] x      The coordinates of the atoms
  */
 void put_atoms_in_box(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x);
 
+/*! \brief Parallellizes put_atoms_in_box()
+ *
+ * This wrapper function around put_atoms_in_box() with the ugly manual
+ * workload splitting is needed to avoid silently introducing multithreading
+ * in tools.
+ * \param[in]     ePBC       The pbc type
+ * \param[in]     box        The simulation box
+ * \param[in,out] x          The coordinates of the atoms
+ * \param[in]     nth        number of threads to be used in the given module
+ */
+void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
+
 /*! \brief Put atoms inside triclinic box
  *
  * This puts ALL atoms in the triclinic unit cell, centered around the
index 7673365bc8fa284335d4147052ab6d05445e2305..68302f8b8f7db31c4d2b6e1678d110fbd3ebba59 100644 (file)
@@ -60,7 +60,6 @@
 #include "gromacs/gmxlib/network.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/mdlib/groupcoord.h"
-#include "gromacs/mdlib/sim_util.h"
 #include "gromacs/mdtypes/commrec.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"