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) 2011,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.
37 #ifndef GMX_TOPOLOGY_TOPOLOGY_H
38 #define GMX_TOPOLOGY_TOPOLOGY_H
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/topology/atoms.h"
46 #include "gromacs/topology/block.h"
47 #include "gromacs/topology/forcefieldparameters.h"
48 #include "gromacs/topology/idef.h"
49 #include "gromacs/topology/symtab.h"
50 #include "gromacs/utility/enumerationhelpers.h"
51 #include "gromacs/utility/unique_cptr.h"
53 enum class SimulationAtomGroupType
: int
61 MassCenterVelocityRemoval
,
62 CompressedPositionOutput
,
63 OrientationRestraintsFit
,
68 //! Short strings used for describing atom groups in log and energy files
69 const char *shortName(SimulationAtomGroupType type
);
71 //const char *shortName(int type); // if necessary
73 /*! \brief Molecules type data: atoms, interactions and exclusions */
80 /*! \brief Deleted copy assignment operator to avoid (not) freeing pointers */
81 gmx_moltype_t
&operator=(const gmx_moltype_t
&) = delete;
83 /*! \brief Default copy constructor */
84 gmx_moltype_t(const gmx_moltype_t
&) = default;
86 char **name
; /**< Name of the molecule type */
87 t_atoms atoms
; /**< The atoms in this molecule */
88 InteractionLists ilist
; /**< Interaction list with local indices */
89 t_block cgs
; /**< The charge groups */
90 t_blocka excls
; /**< The exclusions */
93 /*! \brief Block of molecules of the same type, used in gmx_mtop_t */
96 int type
= -1; /**< The molecule type index in mtop.moltype */
97 int nmol
= 0; /**< The number of molecules in this block */
98 std::vector
<gmx::RVec
> posres_xA
; /**< Position restraint coordinates for top A */
99 std::vector
<gmx::RVec
> posres_xB
; /**< Position restraint coordinates for top B */
102 /*! \brief Indices for a gmx_molblock_t, derived from other gmx_mtop_t contents */
103 struct MoleculeBlockIndices
105 int numAtomsPerMolecule
; /**< Number of atoms in a molecule in the block */
106 int globalAtomStart
; /**< Global atom index of the first atom in the block */
107 int globalAtomEnd
; /**< Global atom index + 1 of the last atom in the block */
108 int globalResidueStart
; /**< Global residue index of the first residue in the block */
109 int residueNumberStart
; /**< Residue numbers start from this value if the number of residues per molecule is <= maxres_renum */
110 int moleculeIndexStart
; /**< Global molecule indexing starts from this value */
113 /*! \brief Contains the simulation atom groups.
115 * Organized as containers for the different
116 * SimulationAtomGroupTypes
118 struct SimulationGroups
120 //! Groups of particles
121 gmx::EnumerationArray
<SimulationAtomGroupType
, AtomGroupIndices
> groups
;
122 //! Names of groups, stored as pointer to the entries in the symbol table.
123 std::vector
<char **> groupNames
;
124 //! Group numbers for the different SimulationAtomGroupType groups.
125 gmx::EnumerationArray
< SimulationAtomGroupType
, std::vector
< unsigned char>> groupNumbers
;
128 * Number of group numbers for a single SimulationGroup.
130 * \param[in] group Integer value for the group type.
132 int numberOfGroupNumbers(SimulationAtomGroupType group
) const { return gmx::ssize(groupNumbers
[group
]); }
136 * Returns group number of an input group for a given atom.
138 * Returns the group \p type for \p atom in \p group, or 0 if the
139 * entries for all atoms in the group are 0 and the pointer is thus null.
141 * \param[in] group Group to check.
142 * \param[in] type Type of group to check.
143 * \param[in] atom Atom to check if it has an entry.
145 int getGroupType (const SimulationGroups
&group
, SimulationAtomGroupType type
, int atom
);
147 /* The global, complete system topology struct, based on molecule types.
148 * This structure should contain no data that is O(natoms) in memory.
150 * TODO: Find a solution for ensuring that the derived data is in sync
151 * with the primary data, possibly by converting to a class.
153 struct gmx_mtop_t
//NOLINT(clang-analyzer-optin.performance.Padding)
159 //! Name of the topology.
160 char **name
= nullptr;
161 //! Force field parameters used.
162 gmx_ffparams_t ffparams
;
163 //! Vector of different molecule types.
164 std::vector
<gmx_moltype_t
> moltype
;
165 //! Vector of different molecule blocks.
166 std::vector
<gmx_molblock_t
> molblock
;
167 //! Are there intermolecular interactions?
168 bool bIntermolecularInteractions
= false;
170 * List of intermolecular interactions using system wide
171 * atom indices, either NULL or size F_NRE
173 std::unique_ptr
<InteractionLists
> intermolecular_ilist
= nullptr;
174 //! Number of global atoms.
176 //! Parameter for residue numbering.
177 int maxres_renum
= 0;
178 //! The maximum residue number in moltype
180 //! Atomtype properties
181 t_atomtypes atomtypes
;
182 //! Groups of atoms for different purposes
183 SimulationGroups groups
;
186 //! Tells whether we have valid molecule indices
187 bool haveMoleculeIndices
= false;
188 /*! \brief List of global atom indices of atoms between which
189 * non-bonded interactions must be excluded.
191 std::vector
<int> intermolecularExclusionGroup
;
193 /* Derived data below */
194 //! Indices for each molblock entry for fast lookup of atom properties
195 std::vector
<MoleculeBlockIndices
> moleculeBlockIndices
;
199 * The fully written out topology for a domain over its lifetime
201 * Also used in some analysis code.
203 struct gmx_localtop_t
205 //! Constructor used for normal operation, manages own resources.
210 //! The interaction function definition
212 //! Atomtype properties
213 t_atomtypes atomtypes
;
216 //! Flag for domain decomposition so we don't free already freed memory.
217 bool useInDomainDecomp_
= false;
220 /* The old topology struct, completely written out, used in analysis tools */
221 typedef struct t_topology
223 char **name
; /* Name of the topology */
224 t_idef idef
; /* The interaction function definition */
225 t_atoms atoms
; /* The atoms */
226 t_atomtypes atomtypes
; /* Atomtype properties */
227 t_block cgs
; /* The charge groups */
228 t_block mols
; /* The molecules */
229 gmx_bool bIntermolecularInteractions
; /* Inter.mol. int. ? */
230 t_blocka excls
; /* The exclusions */
231 t_symtab symtab
; /* The symbol table */
234 void init_top(t_topology
*top
);
235 void done_top(t_topology
*top
);
236 // Frees both t_topology and gmx_mtop_t when the former has been created from
238 void done_top_mtop(t_topology
*top
, gmx_mtop_t
*mtop
);
240 bool gmx_mtop_has_masses(const gmx_mtop_t
*mtop
);
241 bool gmx_mtop_has_charges(const gmx_mtop_t
*mtop
);
242 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t
&mtop
);
243 bool gmx_mtop_has_atomtypes(const gmx_mtop_t
*mtop
);
244 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t
*mtop
);
246 void pr_mtop(FILE *fp
, int indent
, const char *title
, const gmx_mtop_t
*mtop
,
247 gmx_bool bShowNumbers
, gmx_bool bShowParameters
);
248 void pr_top(FILE *fp
, int indent
, const char *title
, const t_topology
*top
,
249 gmx_bool bShowNumbers
, gmx_bool bShowParameters
);
251 /*! \brief Compare two mtop topologies.
253 * \param[in] fp File pointer to write to.
254 * \param[in] mtop1 First topology to compare.
255 * \param[in] mtop2 Second topology to compare.
256 * \param[in] relativeTolerance Relative tolerance for comparison.
257 * \param[in] absoluteTolerance Absolute tolerance for comparison.
259 void compareMtop(FILE *fp
, const gmx_mtop_t
&mtop1
, const gmx_mtop_t
&mtop2
, real relativeTolerance
, real absoluteTolerance
);
261 /*! \brief Check perturbation parameters in topology.
263 * \param[in] fp File pointer to write to.
264 * \param[in] mtop1 Topology to check perturbation parameters in.
265 * \param[in] relativeTolerance Relative tolerance for comparison.
266 * \param[in] absoluteTolerance Absolute tolerance for comparison.
268 void compareMtopAB(FILE *fp
, const gmx_mtop_t
&mtop1
, real relativeTolerance
, real absoluteTolerance
);
270 /*! \brief Compare groups.
272 * \param[in] fp File pointer to write to.
273 * \param[in] g0 First group for comparison.
274 * \param[in] g1 Second group for comparison.
275 * \param[in] natoms0 Number of atoms for first group.
276 * \param[in] natoms1 Number of atoms for second group.
278 void compareAtomGroups(FILE *fp
, const SimulationGroups
&g0
, const SimulationGroups
&g1
,
279 int natoms0
, int natoms1
);
281 //! Typedef for gmx_localtop in analysis tools.
282 using ExpandedTopologyPtr
= std::unique_ptr
<gmx_localtop_t
>;
284 void copy_moltype(const gmx_moltype_t
*src
, gmx_moltype_t
*dst
);