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