d5b633ee061e9f7bae21e03ebed2b4d4d858e0f6
[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 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #ifndef GMX_TOPOLOGY_TOPOLOGY_H
39 #define GMX_TOPOLOGY_TOPOLOGY_H
40
41 #include <cstdio>
42
43 #include <vector>
44
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"
54
55 enum class SimulationAtomGroupType : int
56 {
57     TemperatureCoupling,
58     EnergyOutput,
59     Acceleration,
60     Freeze,
61     User1,
62     User2,
63     MassCenterVelocityRemoval,
64     CompressedPositionOutput,
65     OrientationRestraintsFit,
66     QuantumMechanics,
67     Count
68 };
69
70 //! Short strings used for describing atom groups in log and energy files
71 const char* shortName(SimulationAtomGroupType type);
72
73 // const char *shortName(int type); // if necessary
74
75 /*! \brief Molecules type data: atoms, interactions and exclusions */
76 struct gmx_moltype_t
77 {
78     gmx_moltype_t();
79
80     ~gmx_moltype_t();
81
82     /*! \brief Deleted copy assignment operator to avoid (not) freeing pointers */
83     gmx_moltype_t& operator=(const gmx_moltype_t&) = delete;
84
85     /*! \brief Default copy constructor */
86     gmx_moltype_t(const gmx_moltype_t&) = default;
87
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                       */
92 };
93
94 /*! \brief Block of molecules of the same type, used in gmx_mtop_t */
95 struct gmx_molblock_t
96 {
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 */
101 };
102
103 /*! \brief Indices for a gmx_molblock_t, derived from other gmx_mtop_t contents */
104 struct MoleculeBlockIndices
105 {
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 */
112 };
113
114 /*! \brief Contains the simulation atom groups.
115  *
116  * Organized as containers for the different
117  * SimulationAtomGroupTypes
118  */
119 struct SimulationGroups
120 {
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;
127
128     /*! \brief
129      * Number of group numbers for a single SimulationGroup.
130      *
131      * \param[in] group Integer value for the group type.
132      */
133     int numberOfGroupNumbers(SimulationAtomGroupType group) const
134     {
135         return static_cast<int>(groupNumbers[group].size());
136     }
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     gmx::ListOfLists<int> 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     mols;                        /* The molecules                        */
232     gmx_bool    bIntermolecularInteractions; /* Inter.mol. int. ?   */
233     /* Note that the exclusions are not stored in t_topology */
234     t_symtab symtab; /* The symbol table                     */
235 } t_topology;
236
237 void init_top(t_topology* top);
238 void done_top(t_topology* top);
239 // Frees both t_topology and gmx_mtop_t when the former has been created from
240 // the latter.
241 void done_top_mtop(t_topology* top, gmx_mtop_t* mtop);
242
243 bool gmx_mtop_has_masses(const gmx_mtop_t* mtop);
244 bool gmx_mtop_has_charges(const gmx_mtop_t* mtop);
245 bool gmx_mtop_has_perturbed_charges(const gmx_mtop_t& mtop);
246 bool gmx_mtop_has_atomtypes(const gmx_mtop_t* mtop);
247 bool gmx_mtop_has_pdbinfo(const gmx_mtop_t* mtop);
248
249 void pr_mtop(FILE* fp, int indent, const char* title, const gmx_mtop_t* mtop, gmx_bool bShowNumbers, gmx_bool bShowParameters);
250 void pr_top(FILE* fp, int indent, const char* title, const t_topology* top, gmx_bool bShowNumbers, gmx_bool bShowParameters);
251
252 /*! \brief Compare two mtop topologies.
253  *
254  * \param[in] fp File pointer to write to.
255  * \param[in] mtop1 First topology to compare.
256  * \param[in] mtop2 Second topology to compare.
257  * \param[in] relativeTolerance Relative tolerance for comparison.
258  * \param[in] absoluteTolerance Absolute tolerance for comparison.
259  */
260 void compareMtop(FILE* fp, const gmx_mtop_t& mtop1, const gmx_mtop_t& mtop2, real relativeTolerance, real absoluteTolerance);
261
262 /*! \brief Check perturbation parameters in topology.
263  *
264  * \param[in] fp File pointer to write to.
265  * \param[in] mtop1 Topology to check perturbation parameters in.
266  * \param[in] relativeTolerance Relative tolerance for comparison.
267  * \param[in] absoluteTolerance Absolute tolerance for comparison.
268  */
269 void compareMtopAB(FILE* fp, const gmx_mtop_t& mtop1, real relativeTolerance, real absoluteTolerance);
270
271 /*! \brief Compare groups.
272  *
273  * \param[in] fp File pointer to write to.
274  * \param[in] g0 First group for comparison.
275  * \param[in] g1 Second group for comparison.
276  * \param[in] natoms0 Number of atoms for first group.
277  * \param[in] natoms1 Number of atoms for second group.
278  */
279 void compareAtomGroups(FILE* fp, const SimulationGroups& g0, const SimulationGroups& g1, int natoms0, int natoms1);
280
281 //! Typedef for gmx_localtop in analysis tools.
282 using ExpandedTopologyPtr = std::unique_ptr<gmx_localtop_t>;
283
284 void copy_moltype(const gmx_moltype_t* src, gmx_moltype_t* dst);
285
286 #endif