1 /* 2 * This file is part of the GROMACS molecular simulation package. 3 * 4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands. 5 * Copyright (c) 2001-2004, The GROMACS development team. 6 * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team. 7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by 8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, 9 * and including many others, as listed in the AUTHORS file in the 10 * top-level source directory and at http://www.gromacs.org. 11 * 12 * GROMACS is free software; you can redistribute it and/or 13 * modify it under the terms of the GNU Lesser General Public License 14 * as published by the Free Software Foundation; either version 2.1 15 * of the License, or (at your option) any later version. 16 * 17 * GROMACS is distributed in the hope that it will be useful, 18 * but WITHOUT ANY WARRANTY; without even the implied warranty of 19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 20 * Lesser General Public License for more details. 21 * 22 * You should have received a copy of the GNU Lesser General Public 23 * License along with GROMACS; if not, see 24 * http://www.gnu.org/licenses, or write to the Free Software Foundation, 25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. 26 * 27 * If you want to redistribute modifications to GROMACS, please 28 * consider that scientific software is very special. Version 29 * control is crucial - bugs must be traceable. We will be happy to 30 * consider code for inclusion in the official distribution, but 31 * derived work must not be called official GROMACS. Details are found 32 * in the README & COPYING files - if they are missing, get the 33 * official version at http://www.gromacs.org. 34 * 35 * To help us fund GROMACS development, we humbly ask that you cite 36 * the research papers on the package. Check out http://www.gromacs.org. 37 */ 38 /*! \libinternal \file 39 * \brief Declares the VirtualSitesHandler class and vsite standalone functions 40 * 41 * \author Berk Hess <hess@kth.se> 42 * \ingroup module_mdlib 43 * \inlibraryapi 44 */ 45 46 #ifndef GMX_MDLIB_VSITE_H 47 #define GMX_MDLIB_VSITE_H 48 49 #include <memory> 50 51 #include "gromacs/math/vectypes.h" 52 #include "gromacs/topology/idef.h" 53 #include "gromacs/utility/arrayref.h" 54 #include "gromacs/utility/basedefinitions.h" 55 #include "gromacs/utility/classhelpers.h" 56 #include "gromacs/utility/real.h" 57 58 struct gmx_domdec_t; 59 struct gmx_mtop_t; 60 struct t_commrec; 61 struct InteractionList; 62 struct t_mdatoms; 63 struct t_nrnb; 64 struct gmx_wallcycle; 65 enum class PbcType : int; 66 67 namespace gmx 68 { 69 class RangePartitioning; 70 71 /*! \brief The start value of the vsite indices in the ftype enum 72 * 73 * The validity of the start and end values is checked in makeVirtualSitesHandler(). 74 * This is used to avoid loops over all ftypes just to get the vsite entries. 75 * (We should replace the fixed ilist array by only the used entries.) 76 */ 77 static constexpr int c_ftypeVsiteStart = F_VSITE1; 78 //! The start and end value of the vsite indices in the ftype enum 79 static constexpr int c_ftypeVsiteEnd = F_VSITEN + 1; 80 81 //! Type for storing PBC atom information for all vsite types in the system 82 typedef std::array<std::vector<int>, c_ftypeVsiteEnd - c_ftypeVsiteStart> VsitePbc; 83 84 /*! \libinternal 85 * \brief Class that handles construction of vsites and spreading of vsite forces 86 */ 87 class VirtualSitesHandler 88 { 89 public: 90 //! Constructor, used only be the makeVirtualSitesHandler() factory function 91 VirtualSitesHandler(const gmx_mtop_t& mtop, gmx_domdec_t* domdec, PbcType pbcType); 92 93 ~VirtualSitesHandler(); 94 95 //! Returns the number of virtual sites acting over multiple update groups 96 int numInterUpdategroupVirtualSites() const; 97 98 //! Set VSites and distribute VSite work over threads, should be called after each DD partitioning 99 void setVirtualSites(ArrayRef<const InteractionList> ilist, const t_mdatoms& mdatoms); 100 101 /*! \brief Create positions of vsite atoms based for the local system 102 * 103 * \param[in,out] x The coordinates 104 * \param[in] dt The time step 105 * \param[in,out] v When not empty, velocities for vsites are set as displacement/dt 106 * \param[in] box The box 107 */ 108 void construct(ArrayRef<RVec> x, real dt, ArrayRef<RVec> v, const matrix box) const; 109 110 //! Tells how to handle virial contributions due to virtual sites 111 enum class VirialHandling : int 112 { 113 None, //!< Do not compute virial contributions 114 Pbc, //!< Add contributions working over PBC to shift forces 115 NonLinear //!< Compute contributions due to non-linear virtual sites 116 }; 117 118 /*! \brief Spread the force operating on the vsite atoms on the surrounding atoms. 119 * 120 * vsite should point to a valid object. 121 * The virialHandling parameter determines how virial contributions are handled. 122 * If this is set to Linear, shift forces are accumulated into fshift. 123 * If this is set to NonLinear, non-linear contributions are added to virial. 124 * This non-linear correction is required when the virial is not calculated 125 * afterwards from the particle position and forces, but in a different way, 126 * as for instance for the PME mesh contribution. 127 */ 128 void spreadForces(ArrayRef<const RVec> x, 129 ArrayRef<RVec> f, 130 VirialHandling virialHandling, 131 ArrayRef<RVec> fshift, 132 matrix virial, 133 t_nrnb* nrnb, 134 const matrix box, 135 gmx_wallcycle* wcycle); 136 137 private: 138 //! Implementation type. 139 class Impl; 140 //! Implementation object. 141 PrivateImplPointer<Impl> impl_; 142 }; 143 144 /*! \brief Create positions of vsite atoms based for the local system 145 * 146 * \param[in,out] x The coordinates 147 * \param[in] ip Interaction parameters 148 * \param[in] ilist The interaction list 149 */ 150 void constructVirtualSites(ArrayRef<RVec> x, 151 ArrayRef<const t_iparams> ip, 152 ArrayRef<const InteractionList> ilist); 153 154 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex 155 * 156 * \param[in] mtop The global topology 157 * \param[in,out] x The global coordinates 158 */ 159 void constructVirtualSitesGlobal(const gmx_mtop_t& mtop, ArrayRef<RVec> x); 160 161 //! Tells how to handle virial contributions due to virtual sites 162 enum class VirtualSiteVirialHandling : int 163 { 164 None, //!< Do not compute virial contributions 165 Pbc, //!< Add contributions working over PBC to shift forces 166 NonLinear //!< Compute contributions due to non-linear virtual sites 167 }; 168 169 //! Return the number of non-linear virtual site constructions in the system 170 int countNonlinearVsites(const gmx_mtop_t& mtop); 171 172 /*! \brief Return the number of virtual sites that cross update groups 173 * 174 * \param[in] mtop The global topology 175 * \param[in] updateGroupingPerMoleculetype Update grouping per molecule type, pass empty when not using update groups 176 */ 177 int countInterUpdategroupVsites(const gmx_mtop_t& mtop, 178 ArrayRef<const RangePartitioning> updateGroupingPerMoleculetype); 179 180 /*! \brief Create the virtual site handler 181 * 182 * \param[in] mtop The global topology 183 * \param[in] cr The communication record 184 * \param[in] pbcType The type of PBC 185 * \returns A valid vsite handler object or nullptr when there are no virtual sites 186 */ 187 std::unique_ptr<VirtualSitesHandler> makeVirtualSitesHandler(const gmx_mtop_t& mtop, 188 const t_commrec* cr, 189 PbcType pbcType); 190 191 } // namespace gmx 192 193 #endif 194