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