Update instructions in containers.rst
[gromacs.git] / src / gromacs / topology / topology.h
blobcda4ec7e942e5bba57d907e69566725d709c82f1
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, 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 // TODO: collect groups and groupNumbers in a struct for each group type
122 //! Group numbers for each of the different SimulationAtomGroupType groups.
123 gmx::EnumerationArray<SimulationAtomGroupType, AtomGroupIndices> groups;
124 //! Names of groups, stored as pointer to the entries in the symbol table.
125 std::vector<char**> groupNames;
126 //! Indices into groups for each atom for each of the different SimulationAtomGroupType groups.
127 gmx::EnumerationArray<SimulationAtomGroupType, std::vector<unsigned char>> groupNumbers;
129 /*! \brief
130 * Number of atoms for which group numbers are stored for a single SimulationGroup.
132 * \param[in] group The group type.
134 int numberOfGroupNumbers(SimulationAtomGroupType group) const
136 return static_cast<int>(groupNumbers[group].size());
140 /*! \brief
141 * Returns group number of an input group for a given atom.
143 * Returns the group \p type for \p atom in \p group, or 0 if the
144 * entries for all atoms in the group are 0 and the pointer is thus null.
146 * \param[in] group Group to check.
147 * \param[in] type Type of group to check.
148 * \param[in] atom Atom to check if it has an entry.
150 int getGroupType(const SimulationGroups& group, SimulationAtomGroupType type, int atom);
152 /* The global, complete system topology struct, based on molecule types.
153 * This structure should contain no data that is O(natoms) in memory.
155 * TODO: Find a solution for ensuring that the derived data is in sync
156 * with the primary data, possibly by converting to a class.
158 struct gmx_mtop_t //NOLINT(clang-analyzer-optin.performance.Padding)
160 gmx_mtop_t();
162 ~gmx_mtop_t();
164 //! Name of the topology.
165 char** name = nullptr;
166 //! Force field parameters used.
167 gmx_ffparams_t ffparams;
168 //! Vector of different molecule types.
169 std::vector<gmx_moltype_t> moltype;
170 //! Vector of different molecule blocks.
171 std::vector<gmx_molblock_t> molblock;
172 //! Are there intermolecular interactions?
173 bool bIntermolecularInteractions = false;
174 /* \brief
175 * List of intermolecular interactions using system wide
176 * atom indices, either NULL or size F_NRE
178 std::unique_ptr<InteractionLists> intermolecular_ilist = nullptr;
179 //! Number of global atoms.
180 int natoms = 0;
181 //! Atomtype properties
182 t_atomtypes atomtypes;
183 //! Groups of atoms for different purposes
184 SimulationGroups groups;
185 //! The legacy symbol table
186 t_symtab symtab;
187 //! Tells whether we have valid molecule indices
188 bool haveMoleculeIndices = false;
189 /*! \brief List of global atom indices of atoms between which
190 * non-bonded interactions must be excluded.
192 std::vector<int> intermolecularExclusionGroup;
194 //! Maximum number of residues in molecule to trigger renumbering of residues
195 int maxResiduesPerMoleculeToTriggerRenumber() const
197 return maxResiduesPerMoleculeToTriggerRenumber_;
199 //! Maximum residue number that is not renumbered.
200 int maxResNumberNotRenumbered() const { return maxResNumberNotRenumbered_; }
201 /*! \brief Finalize this data structure.
203 * Should be called after generating or reading mtop, to set some compute
204 * intesive variables to avoid N^2 operations later on.
206 * \todo Move into a builder class, once available.
208 void finalize();
210 /* Derived data below */
211 //! Indices for each molblock entry for fast lookup of atom properties
212 std::vector<MoleculeBlockIndices> moleculeBlockIndices;
214 private:
215 //! Build the molblock indices
216 void buildMolblockIndices();
217 //! Maximum number of residues in molecule to trigger renumbering of residues
218 int maxResiduesPerMoleculeToTriggerRenumber_ = 0;
219 //! The maximum residue number in moltype that is not renumbered
220 int maxResNumberNotRenumbered_ = -1;
223 /*! \brief
224 * The fully written out topology for a domain over its lifetime
226 * Also used in some analysis code.
228 struct gmx_localtop_t
230 //! Constructor used for normal operation, manages own resources.
231 gmx_localtop_t(const gmx_ffparams_t& ffparams);
233 //! The interaction function definition
234 InteractionDefinitions idef;
235 //! The exclusions
236 gmx::ListOfLists<int> excls;
239 /* The old topology struct, completely written out, used in analysis tools */
240 typedef struct t_topology
242 char** name; /* Name of the topology */
243 t_idef idef; /* The interaction function definition */
244 t_atoms atoms; /* The atoms */
245 t_atomtypes atomtypes; /* Atomtype properties */
246 t_block mols; /* The molecules */
247 gmx_bool bIntermolecularInteractions; /* Inter.mol. int. ? */
248 /* Note that the exclusions are not stored in t_topology */
249 t_symtab symtab; /* The symbol table */
250 } t_topology;
252 void init_top(t_topology* top);
253 void done_top(t_topology* top);
254 // Frees both t_topology and gmx_mtop_t when the former has been created from
255 // the latter.
256 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop);
258 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop);
259 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop);
260 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop);
261 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop);
262 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop);
264 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters);
265 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters);
267 /*! \brief Compare two mtop topologies.
269 * \param[in] fp File pointer to write to.
270 * \param[in] mtop1 First topology to compare.
271 * \param[in] mtop2 Second topology to compare.
272 * \param[in] relativeTolerance Relative tolerance for comparison.
273 * \param[in] absoluteTolerance Absolute tolerance for comparison.
275 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance);
277 /*! \brief Check perturbation parameters in topology.
279 * \param[in] fp File pointer to write to.
280 * \param[in] mtop1 Topology to check perturbation parameters in.
281 * \param[in] relativeTolerance Relative tolerance for comparison.
282 * \param[in] absoluteTolerance Absolute tolerance for comparison.
284 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance);
286 /*! \brief Compare groups.
288 * \param[in] fp File pointer to write to.
289 * \param[in] g0 First group for comparison.
290 * \param[in] g1 Second group for comparison.
291 * \param[in] natoms0 Number of atoms for first group.
292 * \param[in] natoms1 Number of atoms for second group.
294 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1);
296 //! Typedef for gmx_localtop in analysis tools.
297 using ExpandedTopologyPtr = std::unique_ptr<gmx_localtop_t>;
299 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst);
301 #endif