2451283b1f232a1e3ef9eb7551872e0c47508546
[alexxy/gromacs.git] / src / gromacs / math / coordinatetransformation.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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 /*! \libinternal \file
36  * \brief Declares coordinate transformation routines.
37  *
38  * \author Christian Blau <blau@kth.se>
39  *
40  * \inlibraryapi
41  * \ingroup module_math
42  */
43 #ifndef GMX_MATH_COORDINATETRANSFORMATION_H
44 #define GMX_MATH_COORDINATETRANSFORMATION_H
45
46 #include "gromacs/math/vectypes.h"
47 #include "gromacs/utility/classhelpers.h"
48
49 #include "matrix.h"
50
51 namespace gmx
52 {
53
54 template<typename>
55 class ArrayRef;
56 class ScaleCoordinates
57 {
58 public:
59     //! Set up coordinate scaling with the scaling factor in each dimension
60     explicit ScaleCoordinates(const RVec& scale);
61     ~ScaleCoordinates();
62
63     //! Copy constructor
64     ScaleCoordinates(const ScaleCoordinates& other);
65     //! Copy assignment
66     ScaleCoordinates& operator=(const ScaleCoordinates& other);
67     //! Move constructor
68     ScaleCoordinates(ScaleCoordinates&& other) noexcept;
69     //! Move assignment
70     ScaleCoordinates& operator=(ScaleCoordinates&& other) noexcept;
71
72     /*! \brief Perform a coordinate transformation on input coordinates.
73      * \param[in] coordinates to be transformed
74      */
75     void operator()(ArrayRef<RVec> coordinates) const;
76
77     /*! \brief Perform a coordinate transformation on an input coordinate.
78      * \param[in] coordinate to be transformed
79      */
80     void operator()(RVec* coordinate) const;
81
82     /*! \brief Apply the inverse scale to coordinates, ignoring dimensions for which scale is zero.
83      * \param[in] coordinates to be transformed
84      */
85     void inverseIgnoringZeroScale(ArrayRef<RVec> coordinates) const;
86
87     /*! \brief Apply the inverse scale to a coordinate, ignoring dimensions for which scale is zero.
88      * \param[in] coordinate to be transformed
89      */
90     void inverseIgnoringZeroScale(RVec* coordinate) const;
91
92 private:
93     class Impl;
94     PrivateImplPointer<Impl> impl_;
95 };
96
97 /*! \libinternal \brief Transform coordinates in three dimensions by first
98  * translating, then scaling them.
99  */
100 class TranslateAndScale
101 {
102 public:
103     /*! \brief Construct a three-dimensional coordinate transformation.
104      * Coordinates are first translated, then scaled.
105      * \param[in] translation to be performed on the coordinates
106      * \param[in] scale to be applied to the coordinates
107      */
108     TranslateAndScale(const RVec& scale, const RVec& translation);
109
110     ~TranslateAndScale();
111
112     //! Copy constructor
113     TranslateAndScale(const TranslateAndScale& other);
114     //! Copy assignment
115     TranslateAndScale& operator=(const TranslateAndScale& other);
116     //! Move constructor
117     TranslateAndScale(TranslateAndScale&& other) noexcept;
118     //! Move assignment
119     TranslateAndScale& operator=(TranslateAndScale&& other) noexcept;
120
121     /*! \brief Perform a coordinate transformation on input coordinates.
122      * \param[in] coordinates to be transformed
123      */
124     void operator()(ArrayRef<RVec> coordinates) const;
125
126     /*! \brief Perform a coordinate transformation on a coordinate.
127      * \param[in] coordinate to be transformed
128      */
129     void operator()(RVec* coordinate) const;
130
131     /*! \brief Returns the scaling operation, discarding the translation.
132      */
133     ScaleCoordinates scaleOperationOnly() const;
134
135 private:
136     class Impl;
137     PrivateImplPointer<Impl> impl_;
138 };
139
140 /*! \libinternal
141  * \brief Affine transformation of three-dimensional coordinates.
142  *
143  * Perfoms in-place coordinate transformations.
144  *
145  * Coordinates are first multiplied by a matrix, then translated.
146  */
147 class AffineTransformation
148 {
149 public:
150     /*! \brief Construct a three-dimensional affine transformation.
151      * \param[in] matrix to be applied to the vectors
152      * \param[in] translation to be performed on the vectors
153      */
154     AffineTransformation(Matrix3x3ConstSpan matrix, const RVec& translation);
155
156     /*! \brief Perform an affine transformation on input vectors.
157      * \param[in,out] vectors to be transformed in-place
158      */
159     void operator()(ArrayRef<RVec> vectors) const;
160
161     /*! \brief Perform an affine transformation on a vector.
162      * \param[in,out] vector to be transformed in-place
163      */
164     void operator()(RVec* vector) const;
165
166 private:
167     //! The matrix describing the affine transformation A(x) = matrix_ * x + translation_
168     Matrix3x3 matrix_;
169     //! The translation vector describing the affine transformation A(x) = matrix * x + translation
170     RVec translation_;
171 };
172
173 } // namespace gmx
174 #endif // CoordinateTransformation