From d854a1395da0af28fc0fceb38611baf28a432818 Mon Sep 17 00:00:00 2001 From: Mark Abraham Date: Mon, 7 Nov 2016 16:58:42 +0100 Subject: [PATCH] Created paddedvector.h This permits forcerec.h to not depend on state.h, and creates a stable location in which to implement e.g. PaddedRVecVector. There are now about 25 fewer source files that depend on state.h. Change-Id: I07c443dea45b394773d60f74b036584af5920cf8 --- src/gromacs/applied-forces/electricfield.cpp | 2 +- src/gromacs/math/CMakeLists.txt | 3 +- src/gromacs/math/paddedvector.h | 65 ++++++++++++++++++++++++++++ src/gromacs/mdtypes/forcerec.h | 7 +-- src/gromacs/mdtypes/state.h | 4 +- 5 files changed, 73 insertions(+), 8 deletions(-) create mode 100644 src/gromacs/math/paddedvector.h diff --git a/src/gromacs/applied-forces/electricfield.cpp b/src/gromacs/applied-forces/electricfield.cpp index 7f6d1d54ca..71f8421f68 100644 --- a/src/gromacs/applied-forces/electricfield.cpp +++ b/src/gromacs/applied-forces/electricfield.cpp @@ -179,7 +179,7 @@ class ElectricField : public IInputRecExtension, public IForceProvider virtual void finishOutput(); virtual void initForcerec(t_forcerec *fr); - // From IForceProvider + //! \copydoc gmx::IForceProvider::calculateForces virtual void calculateForces(const t_commrec *cr, const t_mdatoms *atoms, PaddedRVecVector *force, diff --git a/src/gromacs/math/CMakeLists.txt b/src/gromacs/math/CMakeLists.txt index 181e149961..e57c053241 100644 --- a/src/gromacs/math/CMakeLists.txt +++ b/src/gromacs/math/CMakeLists.txt @@ -1,7 +1,7 @@ # # This file is part of the GROMACS molecular simulation package. # -# Copyright (c) 2013,2014,2015, by the GROMACS development team, led by +# Copyright (c) 2013,2014,2015,2016, 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. @@ -45,6 +45,7 @@ gmx_install_headers( utilities.h vec.h vectypes.h + paddedvector.h ) if (BUILD_TESTING) diff --git a/src/gromacs/math/paddedvector.h b/src/gromacs/math/paddedvector.h new file mode 100644 index 0000000000..7b73dee781 --- /dev/null +++ b/src/gromacs/math/paddedvector.h @@ -0,0 +1,65 @@ +/* + * This file is part of the GROMACS molecular simulation package. + * + * Copyright (c) 2016, 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. + */ +/*! \file + * \brief + * Declares gmx::PaddedRVecVector + * + * \author Mark Abraham + * \inpublicapi + * \ingroup module_math + */ +#ifndef GMX_MATH_PADDEDVECTOR_H +#define GMX_MATH_PADDEDVECTOR_H + +#include + +#include "gromacs/math/vectypes.h" + +namespace gmx +{ + +/*! \brief Temporary definition of a type usable for SIMD-style loads of RVec quantities. + * + * \todo This vector is not padded yet, padding will be added soon */ +using PaddedRVecVector = std::vector; + +} // namespace gmx + +// TODO This is a hack to avoid littering gmx:: all over code that is +// almost all destined to move into the gmx namespace at some point. +// An alternative would be about 20 files with using statements. +using gmx::PaddedRVecVector; + +#endif diff --git a/src/gromacs/mdtypes/forcerec.h b/src/gromacs/mdtypes/forcerec.h index 4b1c0de83d..7cabf08377 100644 --- a/src/gromacs/mdtypes/forcerec.h +++ b/src/gromacs/mdtypes/forcerec.h @@ -38,16 +38,17 @@ #define GMX_MDTYPES_TYPES_FORCEREC_H #include "gromacs/math/vectypes.h" -#include "gromacs/mdtypes/interaction_const.h" -#include "gromacs/mdtypes/md_enums.h" #ifdef __cplusplus -#include "gromacs/mdtypes/state.h" +#include "gromacs/math/paddedvector.h" #endif +#include "gromacs/mdtypes/interaction_const.h" +#include "gromacs/mdtypes/md_enums.h" #include "gromacs/topology/idef.h" #include "gromacs/utility/basedefinitions.h" #include "gromacs/utility/real.h" #ifdef __cplusplus + extern "C" { struct t_mdatoms; diff --git a/src/gromacs/mdtypes/state.h b/src/gromacs/mdtypes/state.h index 475b53e2ca..5b1ed24627 100644 --- a/src/gromacs/mdtypes/state.h +++ b/src/gromacs/mdtypes/state.h @@ -39,6 +39,7 @@ #include +#include "gromacs/math/paddedvector.h" #include "gromacs/math/vectypes.h" #include "gromacs/mdtypes/md_enums.h" #include "gromacs/utility/basedefinitions.h" @@ -71,9 +72,6 @@ enum { /* The names of the state entries, defined in src/gmxlib/checkpoint.c */ extern const char *est_names[estNR]; -/* This vector is not padded yet, padding will be added soon */ -typedef std::vector PaddedRVecVector; - typedef struct history_t { real disre_initf; /* The scaling factor for initializing the time av. */ -- 2.11.4.GIT