2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 /*! \libinternal \file
39 * This file contains the definition of a container for force and virial
42 * Currently the only container defined here is one used in algorithms
43 * that provide their own virial tensor contribution.
44 * We can consider adding another containter for forces and shift forces.
49 * \ingroup module_mdtypes
52 #ifndef GMX_MDTYPES_FORCEOUTPUT_H
53 #define GMX_MDTYPES_FORCEOUTPUT_H
55 #include "gromacs/math/arrayrefwithpadding.h"
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/utility/arrayref.h"
62 /*! \libinternal \brief Container for force and virial for algorithms that compute shift forces for virial calculation
64 * This force output class should be used when computing forces whos virial contribution
65 * is computed using the "single sum virial" algorithm (see the reference manual for
66 * details). To handle the virial contributions of forces working between periodic
67 * images correctly, so-called "shift forces" need to be accumulated for the different
70 class ForceWithShiftForces
73 /*! \brief Constructor
75 * \param[in] force A force buffer that will be used for storing forces
76 * \param[in] computeVirial True when algorithms are required to provide their virial contribution (for the current force evaluation)
77 * \param[in] shiftForces A shift forces buffer of size c_numShiftVectors, only used with \p computeVirial = true
79 ForceWithShiftForces(const gmx::ArrayRefWithPadding<gmx::RVec>& force,
80 const bool computeVirial,
81 const gmx::ArrayRef<gmx::RVec>& shiftForces) :
83 computeVirial_(computeVirial),
84 shiftForces_(computeVirial ? shiftForces : gmx::ArrayRef<gmx::RVec>())
86 GMX_ASSERT(!computeVirial || !shiftForces.empty(),
87 "We need a valid shift force buffer when computing the virial");
90 //! Returns an arrayref to the force buffer without padding
91 gmx::ArrayRef<gmx::RVec> force() { return force_.unpaddedArrayRef(); }
93 //! Returns a const arrayref to the force buffer without padding
94 gmx::ArrayRef<const gmx::RVec> force() const { return force_.unpaddedConstArrayRef(); }
96 //! Returns whether the virial needs to be computed
97 bool computeVirial() const { return computeVirial_; }
99 //! Returns the shift forces buffer
100 gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
102 //! Returns a const shift forces buffer
103 gmx::ArrayRef<const gmx::RVec> shiftForces() const { return shiftForces_; }
105 //! Returns a reference to the boolean which tells whether we have spread forces on vsites
106 bool& haveSpreadVsiteForces() { return haveSpreadVsiteForces_; }
110 gmx::ArrayRefWithPadding<gmx::RVec> force_;
111 //! True when virial computation is requested
113 //! A buffer for storing the shift forces, size c_numShiftVectors
114 gmx::ArrayRef<gmx::RVec> shiftForces_;
115 //! Tells whether we have spread the vsite forces
116 bool haveSpreadVsiteForces_ = false;
119 /*! \libinternal \brief Container for force and virial for algorithms that provide their own virial tensor contribution
121 * \note The \p force_ data member is a reference to an external force buffer.
123 class ForceWithVirial
126 /*! \brief Constructor
128 * \param[in] force A force buffer that will be used for storing forces
129 * \param[in] computeVirial True when algorithms are required to provide their virial contribution (for the current force evaluation)
131 ForceWithVirial(const ArrayRef<RVec>& force, const bool computeVirial) :
132 force_(force), computeVirial_(computeVirial)
134 for (int dim1 = 0; dim1 < DIM; dim1++)
136 for (int dim2 = 0; dim2 < DIM; dim2++)
138 virial_[dim1][dim2] = 0;
143 /*! \brief Adds a virial contribution
145 * \note Can be called with \p computeVirial=false.
146 * \note It is recommended to accumulate the virial contributions
147 * of a module internally before calling this method, as that
148 * will reduce rounding errors.
150 * \param[in] virial The virial contribution to add
152 void addVirialContribution(const matrix virial)
156 for (int dim1 = 0; dim1 < DIM; dim1++)
158 for (int dim2 = 0; dim2 < DIM; dim2++)
160 virial_[dim1][dim2] += virial[dim1][dim2];
166 /*! \brief Adds a virial diagonal contribution
168 * \note Can be called with \p computeVirial=false.
169 * \note It is recommended to accumulate the virial contributions
170 * of a module internally before calling this method, as that
171 * will reduce rounding errors.
173 * \param[in] virial The virial contribution to add
175 void addVirialContribution(const RVec virial)
179 for (int dim = 0; dim < DIM; dim++)
181 virial_[dim][dim] += virial[dim];
186 /*! \brief Returns the accumulated virial contributions
188 const matrix& getVirial() const { return virial_; }
190 const ArrayRef<RVec> force_; //!< Force accumulation buffer reference
191 const bool computeVirial_; //!< True when algorithms are required to provide their virial contribution (for the current force evaluation)
193 matrix virial_; //!< Virial accumulation buffer
196 /*! \libinternal \brief Force and virial output buffers for use in force computation
202 ForceOutputs(const ForceWithShiftForces& forceWithShiftForces,
203 bool haveForceWithVirial,
204 const ForceWithVirial& forceWithVirial) :
205 forceWithShiftForces_(forceWithShiftForces),
206 haveForceWithVirial_(haveForceWithVirial),
207 forceWithVirial_(forceWithVirial)
211 //! Returns a reference to the force with shift forces object
212 ForceWithShiftForces& forceWithShiftForces() { return forceWithShiftForces_; }
214 //! Return whether there are forces with direct virial contributions
215 bool haveForceWithVirial() const { return haveForceWithVirial_; }
217 //! Returns a reference to the force with virial object
218 ForceWithVirial& forceWithVirial() { return forceWithVirial_; }
221 //! Force output buffer used by legacy modules (without SIMD padding)
222 ForceWithShiftForces forceWithShiftForces_;
223 //! Whether we have forces with direct virial contributions
224 bool haveForceWithVirial_;
225 //! Force with direct virial contribution (if there are any; without SIMD padding)
226 ForceWithVirial forceWithVirial_;