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,2019, 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.
38 #ifndef GMX_GMXPREPROCESS_GROMPP_IMPL_H
39 #define GMX_GMXPREPROCESS_GROMPP_IMPL_H
43 #include "gromacs/gmxpreprocess/notset.h"
44 #include "gromacs/topology/atoms.h"
45 #include "gromacs/topology/block.h"
46 #include "gromacs/topology/idef.h"
47 #include "gromacs/utility/arrayref.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 /*! \libinternal \brief
54 * Describes an interaction of a given type, plus its parameters.
56 class InteractionOfType
59 //! Constructor that initializes vectors.
60 InteractionOfType(gmx::ArrayRef<const int> atoms,
61 gmx::ArrayRef<const real> params,
62 const std::string& name = "");
64 //! Access the individual elements set for the parameter.
65 const int& ai() const;
66 const int& aj() const;
67 const int& ak() const;
68 const int& al() const;
69 const int& am() const;
71 const real& c0() const;
72 const real& c1() const;
73 const real& c2() const;
75 const std::string& interactionTypeName() const;
78 /*! \brief Renumbers atom Ids.
80 * Enforces that ai() is less than the opposite terminal atom index,
81 * with the number depending on the interaction type.
85 //! Set single force field parameter.
86 void setForceParameter(int pos, real value);
88 //! View on all atoms numbers that are actually set.
89 gmx::ArrayRef<int> atoms() { return atoms_; }
90 //! Const view on all atoms numbers that are actually set.
91 gmx::ArrayRef<const int> atoms() const { return atoms_; }
92 //! View on all of the force field parameters
93 gmx::ArrayRef<const real> forceParam() const { return forceParam_; }
94 //! View on all of the force field parameters
95 gmx::ArrayRef<real> forceParam() { return forceParam_; }
98 //! Return if we have a bond parameter, means two atoms right now.
99 bool isBond() const { return atoms_.size() == 2; }
100 //! Return if we have an angle parameter, means three atoms right now.
101 bool isAngle() const { return atoms_.size() == 3; }
102 //! Return if we have a dihedral parameter, means four atoms right now.
103 bool isDihedral() const { return atoms_.size() == 4; }
104 //! Return if we have a cmap parameter, means five atoms right now.
105 bool isCmap() const { return atoms_.size() == 5; }
106 //! Enforce that atom id ai() is less than aj().
107 void sortBondAtomIds();
108 //! Enforce that atom id ai() is less than ak(). Does not change aj().
109 void sortAngleAtomIds();
110 /*! \brief Enforce order of atoms in dihedral.
112 * Changes atom order if needed to enforce that ai() is less than al().
113 * If ai() and al() are swapped, aj() and ak() are swapped as well,
114 * independent of their previous order.
116 void sortDihedralAtomIds();
117 //! The atom list (eg. bonds: particle, i = atoms[0] (ai), j = atoms[1] (aj))
118 std::vector<int> atoms_;
119 //! Force parameters (eg. b0 = forceParam[0])
120 std::array<real, MAXFORCEPARAM> forceParam_;
121 //! Used with forcefields whose .rtp files name the interaction types (e.g. GROMOS), rather than look them up from the atom names.
122 std::string interactionTypeName_;
125 /*! \libinternal \brief
126 * A set of interactions of a given type
127 * (found in the enumeration in ifunc.h), complete with
128 * atom indices and force field function parameters.
130 * This is used for containing the data obtained from the
131 * lists of interactions of a given type in a [moleculetype]
132 * topology file definition.
134 struct InteractionsOfType
135 { // NOLINT (clang-analyzer-optin.performance.Padding)
136 //! The different parameters in the system.
137 std::vector<InteractionOfType> interactionTypes;
138 //! CMAP grid spacing.
139 int cmakeGridSpacing = -1;
140 //! Number of cmap angles.
143 std::vector<real> cmap;
144 //! The five atomtypes followed by a number that identifies the type.
145 std::vector<int> cmapAtomTypes;
147 //! Number of parameters.
148 size_t size() const { return interactionTypes.size(); }
149 //! Elements in cmap grid data.
150 int ncmap() const { return cmap.size(); }
151 //! Number of elements in cmapAtomTypes.
152 int nct() const { return cmapAtomTypes.size(); }
157 int nr; /* The number of exclusions */
158 int* e; /* The excluded atoms */
162 /*! \libinternal \brief
163 * Holds the molecule information during preprocessing.
165 struct MoleculeInformation
167 //! Name of the molecule.
168 char** name = nullptr;
169 //! Number of exclusions per atom.
171 //! Have exclusions been generated?.
172 bool excl_set = false;
173 //! Has the mol been processed.
174 bool bProcessed = false;
175 //! Atoms in the moelcule.
177 //! Molecules separated in datastructure.
179 //! Exclusions in the molecule.
180 gmx::ListOfLists<int> excls;
181 //! Interactions of a defined type.
182 std::array<InteractionsOfType, F_NRE> interactions;
187 * This should be removed as soon as the underlying datastructures
188 * have been cleaned up to use proper initialization and can be copy
194 * Partial clean up function.
196 * Should be removed once this datastructure actually owns all its own memory and
197 * elements of it are not stolen by other structures and properly copy constructed
199 * Cleans up the mols and plist datastructures but not cgs and excls.
201 void partialCleanUp();
204 * Full clean up function.
206 * Should be removed once the destructor can always do this.
217 bool is_int(double x);
218 /* Returns TRUE when x is integer */