Use ListOfLists in gmx_mtop_t and gmx_localtop_t
[alexxy/gromacs.git] / src / gromacs / topology / topology.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_TOPOLOGY_TOPOLOGY_H
38 #define GMX_TOPOLOGY_TOPOLOGY_H
39
40 #include <cstdio>
41
42 #include <vector>
43
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/listoflists.h"
52 #include "gromacs/utility/unique_cptr.h"
53
54 enum class SimulationAtomGroupType : int
55 {
56     TemperatureCoupling,
57     EnergyOutput,
58     Acceleration,
59     Freeze,
60     User1,
61     User2,
62     MassCenterVelocityRemoval,
63     CompressedPositionOutput,
64     OrientationRestraintsFit,
65     QuantumMechanics,
66     Count
67 };
68
69 //! Short strings used for describing atom groups in log and energy files
70 const char* shortName(SimulationAtomGroupType type);
71
72 // const char *shortName(int type); // if necessary
73
74 /*! \brief Molecules type data: atoms, interactions and exclusions */
75 struct gmx_moltype_t
76 {
77     gmx_moltype_t();
78
79     ~gmx_moltype_t();
80
81     /*! \brief Deleted copy assignment operator to avoid (not) freeing pointers */
82     gmx_moltype_t& operator=(const gmx_moltype_t&) = delete;
83
84     /*! \brief Default copy constructor */
85     gmx_moltype_t(const gmx_moltype_t&) = default;
86
87     char**                name;  /**< Name of the molecule type            */
88     t_atoms               atoms; /**< The atoms in this molecule           */
89     InteractionLists      ilist; /**< Interaction list with local indices  */
90     gmx::ListOfLists<int> excls; /**< The exclusions                       */
91 };
92
93 /*! \brief Block of molecules of the same type, used in gmx_mtop_t */
94 struct gmx_molblock_t
95 {
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 */
100 };
101
102 /*! \brief Indices for a gmx_molblock_t, derived from other gmx_mtop_t contents */
103 struct MoleculeBlockIndices
104 {
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 */
111 };
112
113 /*! \brief Contains the simulation atom groups.
114  *
115  * Organized as containers for the different
116  * SimulationAtomGroupTypes
117  */
118 struct SimulationGroups
119 {
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;
126
127     /*! \brief
128      * Number of group numbers for a single SimulationGroup.
129      *
130      * \param[in] group Integer value for the group type.
131      */
132     int numberOfGroupNumbers(SimulationAtomGroupType group) const
133     {
134         return gmx::ssize(groupNumbers[group]);
135     }
136 };
137
138 /*! \brief
139  * Returns group number of an input group for a given atom.
140  *
141  * Returns the group \p type for \p atom in \p group, or 0 if the
142  * entries for all atoms in the group are 0 and the pointer is thus null.
143  *
144  * \param[in] group Group to check.
145  * \param[in] type  Type of group to check.
146  * \param[in] atom  Atom to check if it has an entry.
147  */
148 int getGroupType(const SimulationGroups& group, SimulationAtomGroupType type, int atom);
149
150 /* The global, complete system topology struct, based on molecule types.
151  * This structure should contain no data that is O(natoms) in memory.
152  *
153  * TODO: Find a solution for ensuring that the derived data is in sync
154  *       with the primary data, possibly by converting to a class.
155  */
156 struct gmx_mtop_t //NOLINT(clang-analyzer-optin.performance.Padding)
157 {
158     gmx_mtop_t();
159
160     ~gmx_mtop_t();
161
162     //! Name of the topology.
163     char** name = nullptr;
164     //! Force field parameters used.
165     gmx_ffparams_t ffparams;
166     //! Vector of different molecule types.
167     std::vector<gmx_moltype_t> moltype;
168     //! Vector of different molecule blocks.
169     std::vector<gmx_molblock_t> molblock;
170     //! Are there intermolecular interactions?
171     bool bIntermolecularInteractions = false;
172     /* \brief
173      * List of intermolecular interactions using system wide
174      * atom indices, either NULL or size F_NRE
175      */
176     std::unique_ptr<InteractionLists> intermolecular_ilist = nullptr;
177     //! Number of global atoms.
178     int natoms = 0;
179     //! Parameter for residue numbering.
180     int maxres_renum = 0;
181     //! The maximum residue number in moltype
182     int maxresnr = -1;
183     //! Atomtype properties
184     t_atomtypes atomtypes;
185     //! Groups of atoms for different purposes
186     SimulationGroups groups;
187     //! The symbol table
188     t_symtab symtab;
189     //! Tells whether we have valid molecule indices
190     bool haveMoleculeIndices = false;
191     /*! \brief List of global atom indices of atoms between which
192      * non-bonded interactions must be excluded.
193      */
194     std::vector<int> intermolecularExclusionGroup;
195
196     /* Derived data  below */
197     //! Indices for each molblock entry for fast lookup of atom properties
198     std::vector<MoleculeBlockIndices> moleculeBlockIndices;
199 };
200
201 /*! \brief
202  * The fully written out topology for a domain over its lifetime
203  *
204  * Also used in some analysis code.
205  */
206 struct gmx_localtop_t
207 {
208     //! Constructor used for normal operation, manages own resources.
209     gmx_localtop_t();
210
211     ~gmx_localtop_t();
212
213     //! The interaction function definition
214     t_idef idef;
215     //! Atomtype properties
216     t_atomtypes atomtypes;
217     //! The exclusions
218     gmx::ListOfLists<int> excls;
219     //! Flag for domain decomposition so we don't free already freed memory.
220     bool useInDomainDecomp_ = false;
221 };
222
223 /* The old topology struct, completely written out, used in analysis tools */
224 typedef struct t_topology
225 {
226     char**      name;                        /* Name of the topology                 */
227     t_idef      idef;                        /* The interaction function definition  */
228     t_atoms     atoms;                       /* The atoms                            */
229     t_atomtypes atomtypes;                   /* Atomtype properties                  */
230     t_block     mols;                        /* The molecules                        */
231     gmx_bool    bIntermolecularInteractions; /* Inter.mol. int. ?   */
232     t_blocka    excls;                       /* The exclusions                       */
233     t_symtab    symtab;                      /* The symbol table                     */
234 } t_topology;
235
236 void init_top(t_topology* top);
237 void done_top(t_topology* top);
238 // Frees both t_topology and gmx_mtop_t when the former has been created from
239 // the latter.
240 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop);
241
242 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop);
243 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop);
244 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop);
245 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop);
246 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop);
247
248 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters);
249 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters);
250
251 /*! \brief Compare two mtop topologies.
252  *
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.
258  */
259 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance);
260
261 /*! \brief Check perturbation parameters in topology.
262  *
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.
267  */
268 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance);
269
270 /*! \brief Compare groups.
271  *
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.
277  */
278 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1);
279
280 //! Typedef for gmx_localtop in analysis tools.
281 using ExpandedTopologyPtr = std::unique_ptr<gmx_localtop_t>;
282
283 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst);
284
285 #endif