Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdtypes / forceoutput.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \libinternal \file
37  *
38  * \brief
39  * This file contains the definition of a container for force and virial
40  * output.
41  *
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.
45  *
46  * \author Berk Hess
47  *
48  * \inlibraryapi
49  * \ingroup module_mdtypes
50  */
51
52 #ifndef GMX_MDTYPES_FORCEOUTPUT_H
53 #define GMX_MDTYPES_FORCEOUTPUT_H
54
55 #include "gromacs/math/arrayrefwithpadding.h"
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/utility/arrayref.h"
58
59 namespace gmx
60 {
61
62 /*! \libinternal \brief Container for force and virial for algorithms that compute shift forces for virial calculation
63  *
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
68  * periodic images.
69  */
70 class ForceWithShiftForces
71 {
72 public:
73     /*! \brief Constructor
74      *
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 SHIFTS, only used with \p computeVirial = true
78      */
79     ForceWithShiftForces(const gmx::ArrayRefWithPadding<gmx::RVec>& force,
80                          const bool                                 computeVirial,
81                          const gmx::ArrayRef<gmx::RVec>&            shiftForces) :
82         force_(force),
83         computeVirial_(computeVirial),
84         shiftForces_(computeVirial ? shiftForces : gmx::ArrayRef<gmx::RVec>())
85     {
86     }
87
88     //! Returns an arrayref to the force buffer without padding
89     gmx::ArrayRef<gmx::RVec> force() { return force_.unpaddedArrayRef(); }
90
91     //! Returns a const arrayref to the force buffer without padding
92     gmx::ArrayRef<const gmx::RVec> force() const { return force_.unpaddedConstArrayRef(); }
93
94     //! Returns whether the virial needs to be computed
95     bool computeVirial() const { return computeVirial_; }
96
97     //! Returns the shift forces buffer
98     gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
99
100     //! Returns a const shift forces buffer
101     gmx::ArrayRef<const gmx::RVec> shiftForces() const { return shiftForces_; }
102
103 private:
104     //! The force buffer
105     gmx::ArrayRefWithPadding<gmx::RVec> force_;
106     //! True when virial computation is requested
107     bool computeVirial_;
108     //! A buffer for storing the shift forces, size SHIFTS
109     gmx::ArrayRef<gmx::RVec> shiftForces_;
110 };
111
112 /*! \libinternal \brief Container for force and virial for algorithms that provide their own virial tensor contribution
113  *
114  * \note The \p force_ data member is a reference to an external force buffer.
115  */
116 class ForceWithVirial
117 {
118 public:
119     /*! \brief Constructor
120      *
121      * \param[in] force          A force buffer that will be used for storing forces
122      * \param[in] computeVirial  True when algorithms are required to provide their virial contribution (for the current force evaluation)
123      */
124     ForceWithVirial(const ArrayRef<RVec>& force, const bool computeVirial) :
125         force_(force),
126         computeVirial_(computeVirial)
127     {
128         for (int dim1 = 0; dim1 < DIM; dim1++)
129         {
130             for (int dim2 = 0; dim2 < DIM; dim2++)
131             {
132                 virial_[dim1][dim2] = 0;
133             }
134         }
135     }
136
137     /*! \brief Adds a virial contribution
138      *
139      * \note Can be called with \p computeVirial=false.
140      * \note It is recommended to accumulate the virial contributions
141      *       of a module internally before calling this method, as that
142      *       will reduce rounding errors.
143      *
144      * \param[in] virial  The virial contribution to add
145      */
146     void addVirialContribution(const matrix virial)
147     {
148         if (computeVirial_)
149         {
150             for (int dim1 = 0; dim1 < DIM; dim1++)
151             {
152                 for (int dim2 = 0; dim2 < DIM; dim2++)
153                 {
154                     virial_[dim1][dim2] += virial[dim1][dim2];
155                 }
156             }
157         }
158     }
159
160     /*! \brief Adds a virial diagonal contribution
161      *
162      * \note Can be called with \p computeVirial=false.
163      * \note It is recommended to accumulate the virial contributions
164      *       of a module internally before calling this method, as that
165      *       will reduce rounding errors.
166      *
167      * \param[in] virial  The virial contribution to add
168      */
169     void addVirialContribution(const RVec virial)
170     {
171         if (computeVirial_)
172         {
173             for (int dim = 0; dim < DIM; dim++)
174             {
175                 virial_[dim][dim] += virial[dim];
176             }
177         }
178     }
179
180     /*! \brief Returns the accumulated virial contributions
181      */
182     const matrix& getVirial() const { return virial_; }
183
184     const ArrayRef<RVec> force_;         //!< Force accumulation buffer reference
185     const bool           computeVirial_; //!< True when algorithms are required to provide their virial contribution (for the current force evaluation)
186 private:
187     matrix virial_; //!< Virial accumulation buffer
188 };
189
190 /*! \libinternal \brief Force and virial output buffers for use in force computation
191  */
192 class ForceOutputs
193 {
194 public:
195     //! Constructor
196     ForceOutputs(const ForceWithShiftForces& forceWithShiftForces, const ForceWithVirial& forceWithVirial) :
197         forceWithShiftForces_(forceWithShiftForces),
198         forceWithVirial_(forceWithVirial)
199     {
200     }
201
202     //! Returns a reference to the force with shift forces object
203     ForceWithShiftForces& forceWithShiftForces() { return forceWithShiftForces_; }
204
205     //! Returns a reference to the force with virial object
206     ForceWithVirial& forceWithVirial() { return forceWithVirial_; }
207
208 private:
209     //! Force output buffer used by legacy modules (without SIMD padding)
210     ForceWithShiftForces forceWithShiftForces_;
211     //! Force with direct virial contribution (if there are any; without SIMD padding)
212     ForceWithVirial forceWithVirial_;
213 };
214
215 } // namespace gmx
216
217 #endif