2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, 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/vectypes.h"
56 #include "gromacs/utility/arrayref.h"
61 /*! \libinternal \brief Container for force and virial for algorithms that provide their own virial tensor contribution
63 * \note The \p force_ data member is a reference to an external force buffer.
68 /*! \brief Constructor
70 * \param[in] force A force buffer that will be used for storing forces
71 * \param[in] computeVirial True when algorithms are required to provide their virial contribution (for the current force evaluation)
73 ForceWithVirial(ArrayRef
<RVec
> force
, bool computeVirial
) :
75 computeVirial_(computeVirial
)
77 for (int dim1
= 0; dim1
< DIM
; dim1
++)
79 for (int dim2
= 0; dim2
< DIM
; dim2
++)
81 virial_
[dim1
][dim2
] = 0;
86 /*! \brief Adds a virial contribution
88 * \note Can be called with \p computeVirial=false.
89 * \note It is recommended to accumulate the virial contributions
90 * of a module internally before calling this method, as that
91 * will reduce rounding errors.
93 * \param[in] virial The virial contribution to add
95 void addVirialContribution(const matrix virial
)
99 for (int dim1
= 0; dim1
< DIM
; dim1
++)
101 for (int dim2
= 0; dim2
< DIM
; dim2
++)
103 virial_
[dim1
][dim2
] += virial
[dim1
][dim2
];
109 /*! \brief Adds a virial diagonal contribution
111 * \note Can be called with \p computeVirial=false.
112 * \note It is recommended to accumulate the virial contributions
113 * of a module internally before calling this method, as that
114 * will reduce rounding errors.
116 * \param[in] virial The virial contribution to add
118 void addVirialContribution(const RVec virial
)
122 for (int dim
= 0; dim
< DIM
; dim
++)
124 virial_
[dim
][dim
] += virial
[dim
];
129 /*! \brief Returns the accumulated virial contributions
131 const matrix
&getVirial() const
136 const ArrayRef
<RVec
> force_
; //!< Force accumulation buffer reference
137 const bool computeVirial_
; //!< True when algorithms are required to provide their virial contribution (for the current force evaluation)
139 matrix virial_
; //!< Virial accumulation buffer
142 /*! \libinternal \brief Force and virial output buffers for use in force computation
144 * TODO: Extend with shift forces
150 ForceOutputs(rvec
*f
,
151 gmx::ForceWithVirial forceWithVirial
) :
153 forceWithVirial_(forceWithVirial
) {}
155 //! Returns a, deprecated, rvec pointer to the force buffer for use with shift forces
161 //! Returns a reference to the force with virial object
162 gmx::ForceWithVirial
&forceWithVirial()
164 return forceWithVirial_
;
168 //! Force output buffer used by legacy modules (without SIMD padding)
170 //! Force with direct virial contribution (if there are any; without SIMD padding)
171 gmx::ForceWithVirial forceWithVirial_
;