Use moltype directly for obtaining residueAtomRanges
[gromacs.git] / src / gromacs / topology / topology.h
bloba3fe344f77b20e8d0c76de4aef9d053b92e8d0e5
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) 2011,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.
38 #ifndef GMX_TOPOLOGY_TOPOLOGY_H
39 #define GMX_TOPOLOGY_TOPOLOGY_H
41 #include <cstdio>
43 #include <vector>
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/topology/atoms.h"
47 #include "gromacs/topology/block.h"
48 #include "gromacs/topology/forcefieldparameters.h"
49 #include "gromacs/topology/idef.h"
50 #include "gromacs/topology/symtab.h"
51 #include "gromacs/utility/enumerationhelpers.h"
52 #include "gromacs/utility/listoflists.h"
53 #include "gromacs/utility/unique_cptr.h"
55 enum class SimulationAtomGroupType : int
57 TemperatureCoupling,
58 EnergyOutput,
59 Acceleration,
60 Freeze,
61 User1,
62 User2,
63 MassCenterVelocityRemoval,
64 CompressedPositionOutput,
65 OrientationRestraintsFit,
66 QuantumMechanics,
67 Count
70 //! Short strings used for describing atom groups in log and energy files
71 const char* shortName(SimulationAtomGroupType type);
73 // const char *shortName(int type); // if necessary
75 /*! \brief Molecules type data: atoms, interactions and exclusions */
76 struct gmx_moltype_t
78 gmx_moltype_t();
80 ~gmx_moltype_t();
82 /*! \brief Deleted copy assignment operator to avoid (not) freeing pointers */
83 gmx_moltype_t& operator=(const gmx_moltype_t&) = delete;
85 /*! \brief Default copy constructor */
86 gmx_moltype_t(const gmx_moltype_t&) = default;
88 char** name; /**< Name of the molecule type */
89 t_atoms atoms; /**< The atoms in this molecule */
90 InteractionLists ilist; /**< Interaction list with local indices */
91 gmx::ListOfLists<int> excls; /**< The exclusions */
94 /*! \brief Block of molecules of the same type, used in gmx_mtop_t */
95 struct gmx_molblock_t
97 int type = -1; /**< The molecule type index in mtop.moltype */
98 int nmol = 0; /**< The number of molecules in this block */
99 std::vector<gmx::RVec> posres_xA; /**< Position restraint coordinates for top A */
100 std::vector<gmx::RVec> posres_xB; /**< Position restraint coordinates for top B */
103 /*! \brief Indices for a gmx_molblock_t, derived from other gmx_mtop_t contents */
104 struct MoleculeBlockIndices
106 int numAtomsPerMolecule; /**< Number of atoms in a molecule in the block */
107 int globalAtomStart; /**< Global atom index of the first atom in the block */
108 int globalAtomEnd; /**< Global atom index + 1 of the last atom in the block */
109 int globalResidueStart; /**< Global residue index of the first residue in the block */
110 int residueNumberStart; /**< Residue numbers start from this value if the number of residues per molecule is <= maxres_renum */
111 int moleculeIndexStart; /**< Global molecule indexing starts from this value */
114 /*! \brief Contains the simulation atom groups.
116 * Organized as containers for the different
117 * SimulationAtomGroupTypes
119 struct SimulationGroups
121 //! Groups of particles
122 gmx::EnumerationArray<SimulationAtomGroupType, AtomGroupIndices> groups;
123 //! Names of groups, stored as pointer to the entries in the symbol table.
124 std::vector<char**> groupNames;
125 //! Group numbers for the different SimulationAtomGroupType groups.
126 gmx::EnumerationArray<SimulationAtomGroupType, std::vector<unsigned char>> groupNumbers;
128 /*! \brief
129 * Number of group numbers for a single SimulationGroup.
131 * \param[in] group Integer value for the group type.
133 int numberOfGroupNumbers(SimulationAtomGroupType group) const
135 return static_cast<int>(groupNumbers[group].size());
139 /*! \brief
140 * Returns group number of an input group for a given atom.
142 * Returns the group \p type for \p atom in \p group, or 0 if the
143 * entries for all atoms in the group are 0 and the pointer is thus null.
145 * \param[in] group Group to check.
146 * \param[in] type Type of group to check.
147 * \param[in] atom Atom to check if it has an entry.
149 int getGroupType(const SimulationGroups& group, SimulationAtomGroupType type, int atom);
151 /* The global, complete system topology struct, based on molecule types.
152 * This structure should contain no data that is O(natoms) in memory.
154 * TODO: Find a solution for ensuring that the derived data is in sync
155 * with the primary data, possibly by converting to a class.
157 struct gmx_mtop_t //NOLINT(clang-analyzer-optin.performance.Padding)
159 gmx_mtop_t();
161 ~gmx_mtop_t();
163 //! Name of the topology.
164 char** name = nullptr;
165 //! Force field parameters used.
166 gmx_ffparams_t ffparams;
167 //! Vector of different molecule types.
168 std::vector<gmx_moltype_t> moltype;
169 //! Vector of different molecule blocks.
170 std::vector<gmx_molblock_t> molblock;
171 //! Are there intermolecular interactions?
172 bool bIntermolecularInteractions = false;
173 /* \brief
174 * List of intermolecular interactions using system wide
175 * atom indices, either NULL or size F_NRE
177 std::unique_ptr<InteractionLists> intermolecular_ilist = nullptr;
178 //! Number of global atoms.
179 int natoms = 0;
180 //! Parameter for residue numbering.
181 int maxres_renum = 0;
182 //! The maximum residue number in moltype
183 int maxresnr = -1;
184 //! Atomtype properties
185 t_atomtypes atomtypes;
186 //! Groups of atoms for different purposes
187 SimulationGroups groups;
188 //! The legacy symbol table
189 t_symtab symtab;
190 //! Tells whether we have valid molecule indices
191 bool haveMoleculeIndices = false;
192 /*! \brief List of global atom indices of atoms between which
193 * non-bonded interactions must be excluded.
195 std::vector<int> intermolecularExclusionGroup;
197 /* Derived data below */
198 //! Indices for each molblock entry for fast lookup of atom properties
199 std::vector<MoleculeBlockIndices> moleculeBlockIndices;
202 /*! \brief
203 * The fully written out topology for a domain over its lifetime
205 * Also used in some analysis code.
207 struct gmx_localtop_t
209 //! Constructor used for normal operation, manages own resources.
210 gmx_localtop_t(const gmx_ffparams_t& ffparams);
212 //! The interaction function definition
213 InteractionDefinitions idef;
214 //! The exclusions
215 gmx::ListOfLists<int> excls;
218 /* The old topology struct, completely written out, used in analysis tools */
219 typedef struct t_topology
221 char** name; /* Name of the topology */
222 t_idef idef; /* The interaction function definition */
223 t_atoms atoms; /* The atoms */
224 t_atomtypes atomtypes; /* Atomtype properties */
225 t_block mols; /* The molecules */
226 gmx_bool bIntermolecularInteractions; /* Inter.mol. int. ? */
227 /* Note that the exclusions are not stored in t_topology */
228 t_symtab symtab; /* The symbol table */
229 } t_topology;
231 void init_top(t_topology* top);
232 void done_top(t_topology* top);
233 // Frees both t_topology and gmx_mtop_t when the former has been created from
234 // the latter.
235 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop);
237 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop);
238 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop);
239 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop);
240 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop);
241 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop);
243 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters);
244 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters);
246 /*! \brief Compare two mtop topologies.
248 * \param[in] fp File pointer to write to.
249 * \param[in] mtop1 First topology to compare.
250 * \param[in] mtop2 Second topology to compare.
251 * \param[in] relativeTolerance Relative tolerance for comparison.
252 * \param[in] absoluteTolerance Absolute tolerance for comparison.
254 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance);
256 /*! \brief Check perturbation parameters in topology.
258 * \param[in] fp File pointer to write to.
259 * \param[in] mtop1 Topology to check perturbation parameters in.
260 * \param[in] relativeTolerance Relative tolerance for comparison.
261 * \param[in] absoluteTolerance Absolute tolerance for comparison.
263 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance);
265 /*! \brief Compare groups.
267 * \param[in] fp File pointer to write to.
268 * \param[in] g0 First group for comparison.
269 * \param[in] g1 Second group for comparison.
270 * \param[in] natoms0 Number of atoms for first group.
271 * \param[in] natoms1 Number of atoms for second group.
273 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1);
275 //! Typedef for gmx_localtop in analysis tools.
276 using ExpandedTopologyPtr = std::unique_ptr<gmx_localtop_t>;
278 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst);
280 #endif