Parse user-supplied GPU task assignment only when needed
[gromacs.git] / src / gromacs / math / vectypes.h
blob71511750c720a5ea83e1c4404806f5956112e39c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2016,2017, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_MATH_VECTYPES_H
38 #define GMX_MATH_VECTYPES_H
40 #include "gromacs/utility/real.h"
42 #define XX 0 /* Defines for indexing in */
43 #define YY 1 /* vectors */
44 #define ZZ 2
45 #define DIM 3 /* Dimension of vectors */
47 typedef real rvec[DIM];
49 typedef double dvec[DIM];
51 typedef real matrix[DIM][DIM];
53 typedef real tensor[DIM][DIM];
55 typedef int ivec[DIM];
57 #ifdef __cplusplus
59 namespace gmx
62 /*! \brief
63 * C++ class for 3D vectors.
65 * \tparam ValueType Type
67 * This class provides a C++ version of rvec/dvec/ivec that can be put into STL
68 * containers etc. It is more or less a drop-in replacement for `rvec` and
69 * friends: it can be used in most contexts that accept the equivalent C type.
70 * However, there are two cases where explicit conversion is necessary:
71 * - An array of these objects needs to be converted with as_vec_array() (or
72 * convenience methods like as_rvec_array()).
73 * - Passing an RVec as a `const rvec &` parameter to a function needs an
74 * explicit call to as_vec(). The implicit conversion should work for this
75 * as well, but cppcheck parses the necessary implicit conversion operator
76 * incorrectly and MSVC fails to compile code that relies on the implicit
77 * conversion, so the explicit method is necessary.
79 * For the array conversion to work, the compiler should not add any extra
80 * alignment/padding in the layout of this class; that this actually works as
81 * intended is tested in the unit tests.
83 * \inpublicapi
85 template <typename ValueType>
86 class BasicVector
88 public:
89 //! Underlying raw C array type (rvec/dvec/ivec).
90 typedef ValueType RawArray[DIM];
92 //! Constructs default (uninitialized) vector.
93 BasicVector() {}
94 //! Constructs a vector from given values.
95 BasicVector(ValueType x, ValueType y, ValueType z)
97 x_[XX] = x;
98 x_[YY] = y;
99 x_[ZZ] = z;
101 /*! \brief
102 * Constructs a vector from given values.
104 * This constructor is not explicit to support implicit conversions
105 * that allow, e.g., calling `std::vector<RVec>:``:push_back()` directly
106 * with an `rvec` parameter.
108 BasicVector(const RawArray x)
110 x_[XX] = x[XX];
111 x_[YY] = x[YY];
112 x_[ZZ] = x[ZZ];
114 //! Indexing operator to make the class work as the raw array.
115 ValueType &operator[](int i) { return x_[i]; }
116 //! Indexing operator to make the class work as the raw array.
117 ValueType operator[](int i) const { return x_[i]; }
118 // The conversion functions below could more accurately return
119 // RawArray &, but this fails with cppcheck and does not solve the
120 // issue with MSVC, so as_vec() should be used instead.
121 //! Makes BasicVector usable in contexts where a raw C array is expected.
122 operator ValueType *() { return x_; }
123 //! Makes BasicVector usable in contexts where a raw C array is expected.
124 operator const ValueType *() const { return x_; }
126 //! Converts to a raw C array where implicit conversion does not work.
127 RawArray &as_vec() { return x_; }
128 //! Converts to a raw C array where implicit conversion does not work.
129 const RawArray &as_vec() const { return x_; }
131 private:
132 RawArray x_;
135 /*! \brief
136 * Casts a gmx::BasicVector array into an equivalent raw C array.
138 template <typename ValueType> static inline
139 typename BasicVector<ValueType>::RawArray *
140 as_vec_array(BasicVector<ValueType> *x)
142 return reinterpret_cast<typename BasicVector<ValueType>::RawArray *>(x);
145 /*! \brief
146 * Casts a gmx::BasicVector array into an equivalent raw C array.
148 template <typename ValueType> static inline
149 const typename BasicVector<ValueType>::RawArray *
150 as_vec_array(const BasicVector<ValueType> *x)
152 return reinterpret_cast<const typename BasicVector<ValueType>::RawArray *>(x);
155 //! Shorthand for C++ `rvec`-equivalent type.
156 typedef BasicVector<real> RVec;
157 //! Casts a gmx::RVec array into an `rvec` array.
158 static inline rvec *as_rvec_array(RVec *x)
160 return as_vec_array(x);
162 //! Casts a gmx::RVec array into an `rvec` array.
163 static inline const rvec *as_rvec_array(const RVec *x)
165 return as_vec_array(x);
168 //! Shorthand for C++ `ivec`-equivalent type.
169 typedef BasicVector<int> IVec;
171 } // namespace gmx
173 #endif
175 #endif