Remove charge groups from domdec and localtop
[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_block             cgs;    /**< The charge groups                    */
90     t_blocka            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     SimulationGroups();
121
122     ~SimulationGroups();
123
124     //! Groups of particles
125     gmx::EnumerationArray<SimulationAtomGroupType, t_grps>                     groups;
126     //! Names of groups, stored as pointer to the entries in the symbol table.
127     std::vector<char **>                                                       groupNames;
128     //! Group numbers for the different SimulationAtomGroupType groups.
129     gmx::EnumerationArray < SimulationAtomGroupType, std::vector < unsigned char>> groupNumbers;
130
131     /*! \brief
132      * Number of group numbers for a single SimulationGroup.
133      *
134      * \param[in] group Integer value for the group type.
135      */
136     int numberOfGroupNumbers(SimulationAtomGroupType group) const { return groupNumbers[group].size(); }
137 };
138
139 /*! \brief
140  * Returns group number of an input group for a given atom.
141  *
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.
144  *
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.
148  */
149 int getGroupType (const SimulationGroups &group, SimulationAtomGroupType type, int atom);
150
151 /* The global, complete system topology struct, based on molecule types.
152  * This structure should contain no data that is O(natoms) in memory.
153  *
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.
156  */
157 struct gmx_mtop_t //NOLINT(clang-analyzer-optin.performance.Padding)
158 {
159     gmx_mtop_t();
160
161     ~gmx_mtop_t();
162
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
176      */
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 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.
194      */
195     std::vector<int>                  intermolecularExclusionGroup;
196
197     /* Derived data  below */
198     //! Indices for each molblock entry for fast lookup of atom properties
199     std::vector<MoleculeBlockIndices> moleculeBlockIndices;
200 };
201
202 /*! \brief
203  * The fully written out topology for a domain over its lifetime
204  *
205  * Also used in some analysis code.
206  */
207 struct gmx_localtop_t
208 {
209     //! Constructor used for normal operation, manages own resources.
210     gmx_localtop_t();
211
212     ~gmx_localtop_t();
213
214     //! The interaction function definition
215     t_idef        idef;
216     //! Atomtype properties
217     t_atomtypes   atomtypes;
218     //! The exclusions
219     t_blocka      excls;
220     //! Flag for domain decomposition so we don't free already freed memory.
221     bool          useInDomainDecomp_ = false;
222 };
223
224 /* The old topology struct, completely written out, used in analysis tools */
225 typedef struct t_topology
226 {
227     char          **name;                        /* Name of the topology                 */
228     t_idef          idef;                        /* The interaction function definition  */
229     t_atoms         atoms;                       /* The atoms                            */
230     t_atomtypes     atomtypes;                   /* Atomtype properties                  */
231     t_block         cgs;                         /* The charge groups                    */
232     t_block         mols;                        /* The molecules                        */
233     gmx_bool        bIntermolecularInteractions; /* Inter.mol. int. ?   */
234     t_blocka        excls;                       /* The exclusions                       */
235     t_symtab        symtab;                      /* The symbol table                     */
236 } t_topology;
237
238 void init_top(t_topology *top);
239 void done_top(t_topology *top);
240 // Frees both t_topology and gmx_mtop_t when the former has been created from
241 // the latter.
242 void done_top_mtop(t_topology *top, gmx_mtop_t *mtop);
243
244 bool gmx_mtop_has_masses(const gmx_mtop_t *mtop);
245 bool gmx_mtop_has_charges(const gmx_mtop_t *mtop);
246 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t &mtop);
247 bool gmx_mtop_has_atomtypes(const gmx_mtop_t *mtop);
248 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t *mtop);
249
250 void pr_mtop(FILE *fp, int indent, const char *title, const gmx_mtop_t *mtop,
251              gmx_bool bShowNumbers, gmx_bool bShowParameters);
252 void pr_top(FILE *fp, int indent, const char *title, const t_topology *top,
253             gmx_bool bShowNumbers, gmx_bool bShowParameters);
254
255 /*! \brief Compare two mtop topologies.
256  *
257  * \param[in] fp File pointer to write to.
258  * \param[in] mtop1 First topology to compare.
259  * \param[in] mtop2 Second topology to compare.
260  * \param[in] relativeTolerance Relative tolerance for comparison.
261  * \param[in] absoluteTolerance Absolute tolerance for comparison.
262  */
263 void compareMtop(FILE *fp, const gmx_mtop_t &mtop1, const gmx_mtop_t &mtop2, real relativeTolerance, real absoluteTolerance);
264
265 /*! \brief Check perturbation parameters in topology.
266  *
267  * \param[in] fp File pointer to write to.
268  * \param[in] mtop1 Topology to check perturbation parameters in.
269  * \param[in] relativeTolerance Relative tolerance for comparison.
270  * \param[in] absoluteTolerance Absolute tolerance for comparison.
271  */
272 void compareMtopAB(FILE *fp, const gmx_mtop_t &mtop1, real relativeTolerance, real absoluteTolerance);
273
274 /*! \brief Compare groups.
275  *
276  * \param[in] fp File pointer to write to.
277  * \param[in] g0 First group for comparison.
278  * \param[in] g1 Second group for comparison.
279  * \param[in] natoms0 Number of atoms for first group.
280  * \param[in] natoms1 Number of atoms for second group.
281  */
282 void compareAtomGroups(FILE *fp, const SimulationGroups &g0, const SimulationGroups &g1,
283                        int natoms0, int natoms1);
284
285 //! Typedef for gmx_localtop in analysis tools.
286 using ExpandedTopologyPtr = std::unique_ptr<gmx_localtop_t>;
287
288 void copy_moltype(const gmx_moltype_t *src, gmx_moltype_t *dst);
289
290 #endif