Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxpreprocess / grompp_impl.h
blob47ecb013b4c2f99ed0d5a6a796268c7334b695a4
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,2015,2016,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 #ifndef GMX_GMXPREPROCESS_GROMPP_IMPL_H
40 #define GMX_GMXPREPROCESS_GROMPP_IMPL_H
42 #include <string>
44 #include "gromacs/gmxpreprocess/notset.h"
45 #include "gromacs/topology/atoms.h"
46 #include "gromacs/topology/block.h"
47 #include "gromacs/topology/idef.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/listoflists.h"
51 #include "gromacs/utility/real.h"
53 namespace gmx
55 template<typename>
56 class ArrayRef;
59 /*! \libinternal \brief
60 * Describes an interaction of a given type, plus its parameters.
62 class InteractionOfType
64 public:
65 //! Constructor that initializes vectors.
66 InteractionOfType(gmx::ArrayRef<const int> atoms,
67 gmx::ArrayRef<const real> params,
68 const std::string& name = "");
69 /*!@{*/
70 //! Access the individual elements set for the parameter.
71 const int& ai() const;
72 const int& aj() const;
73 const int& ak() const;
74 const int& al() const;
75 const int& am() const;
77 const real& c0() const;
78 const real& c1() const;
79 const real& c2() const;
81 const std::string& interactionTypeName() const;
82 /*!@}*/
84 /*! \brief Renumbers atom Ids.
86 * Enforces that ai() is less than the opposite terminal atom index,
87 * with the number depending on the interaction type.
89 void sortAtomIds();
91 //! Set single force field parameter.
92 void setForceParameter(int pos, real value);
94 //! View on all atoms numbers that are actually set.
95 gmx::ArrayRef<int> atoms() { return atoms_; }
96 //! Const view on all atoms numbers that are actually set.
97 gmx::ArrayRef<const int> atoms() const { return atoms_; }
98 //! View on all of the force field parameters
99 gmx::ArrayRef<const real> forceParam() const { return forceParam_; }
100 //! View on all of the force field parameters
101 gmx::ArrayRef<real> forceParam() { return forceParam_; }
103 private:
104 //! Return if we have a bond parameter, means two atoms right now.
105 bool isBond() const { return atoms_.size() == 2; }
106 //! Return if we have an angle parameter, means three atoms right now.
107 bool isAngle() const { return atoms_.size() == 3; }
108 //! Return if we have a dihedral parameter, means four atoms right now.
109 bool isDihedral() const { return atoms_.size() == 4; }
110 //! Return if we have a cmap parameter, means five atoms right now.
111 bool isCmap() const { return atoms_.size() == 5; }
112 //! Enforce that atom id ai() is less than aj().
113 void sortBondAtomIds();
114 //! Enforce that atom id ai() is less than ak(). Does not change aj().
115 void sortAngleAtomIds();
116 /*! \brief Enforce order of atoms in dihedral.
118 * Changes atom order if needed to enforce that ai() is less than al().
119 * If ai() and al() are swapped, aj() and ak() are swapped as well,
120 * independent of their previous order.
122 void sortDihedralAtomIds();
123 //! The atom list (eg. bonds: particle, i = atoms[0] (ai), j = atoms[1] (aj))
124 std::vector<int> atoms_;
125 //! Force parameters (eg. b0 = forceParam[0])
126 std::array<real, MAXFORCEPARAM> forceParam_;
127 //! Used with forcefields whose .rtp files name the interaction types (e.g. GROMOS), rather than look them up from the atom names.
128 std::string interactionTypeName_;
131 /*! \libinternal \brief
132 * A set of interactions of a given type
133 * (found in the enumeration in ifunc.h), complete with
134 * atom indices and force field function parameters.
136 * This is used for containing the data obtained from the
137 * lists of interactions of a given type in a [moleculetype]
138 * topology file definition.
140 struct InteractionsOfType
141 { // NOLINT (clang-analyzer-optin.performance.Padding)
142 //! The different parameters in the system.
143 std::vector<InteractionOfType> interactionTypes;
144 //! CMAP grid spacing.
145 int cmakeGridSpacing = -1;
146 //! Number of cmap angles.
147 int cmapAngles = -1;
148 //! CMAP grid data.
149 std::vector<real> cmap;
150 //! The five atomtypes followed by a number that identifies the type.
151 std::vector<int> cmapAtomTypes;
153 //! Number of parameters.
154 size_t size() const { return interactionTypes.size(); }
155 //! Elements in cmap grid data.
156 int ncmap() const { return cmap.size(); }
157 //! Number of elements in cmapAtomTypes.
158 int nct() const { return cmapAtomTypes.size(); }
161 struct t_excls
163 int nr; /* The number of exclusions */
164 int* e; /* The excluded atoms */
168 /*! \libinternal \brief
169 * Holds the molecule information during preprocessing.
171 struct MoleculeInformation
173 //! Name of the molecule.
174 char** name = nullptr;
175 //! Number of exclusions per atom.
176 int nrexcl = 0;
177 //! Have exclusions been generated?.
178 bool excl_set = false;
179 //! Has the mol been processed.
180 bool bProcessed = false;
181 //! Atoms in the moelcule.
182 t_atoms atoms;
183 //! Molecules separated in datastructure.
184 t_block mols;
185 //! Exclusions in the molecule.
186 gmx::ListOfLists<int> excls;
187 //! Interactions of a defined type.
188 std::array<InteractionsOfType, F_NRE> interactions;
190 /*! \brief
191 * Initializer.
193 * This should be removed as soon as the underlying datastructures
194 * have been cleaned up to use proper initialization and can be copy
195 * constructed.
197 void initMolInfo();
199 /*! \brief
200 * Partial clean up function.
202 * Should be removed once this datastructure actually owns all its own memory and
203 * elements of it are not stolen by other structures and properly copy constructed
204 * or moved.
205 * Cleans up the mols and plist datastructures but not cgs and excls.
207 void partialCleanUp();
209 /*! \brief
210 * Full clean up function.
212 * Should be removed once the destructor can always do this.
214 void fullCleanUp();
217 struct t_mols
219 std::string name;
220 int nr;
223 bool is_int(double x);
224 /* Returns TRUE when x is integer */
226 #endif