Remove unnecessary includes of arrayref.h
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / hackblock.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) 2013,2014,2015,2018,2019,2020, 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 /*! \file
38  * \libinternal \brief
39  * Methods to modify atoms during preprocessing.
40  */
41 #ifndef GMX_GMXPREPROCESS_HACKBLOCK_H
42 #define GMX_GMXPREPROCESS_HACKBLOCK_H
43
44 #include <cstdio>
45
46 #include <array>
47 #include <string>
48 #include <vector>
49
50 #include "gromacs/gmxpreprocess/notset.h"
51 #include "gromacs/topology/ifunc.h"
52
53 struct t_atom;
54 struct t_symtab;
55
56 namespace gmx
57 {
58 template<typename>
59 class ArrayRef;
60 }
61
62 /*! \brief
63  * Used for reading .rtp/.tdb
64  * ebtsBONDS must be the first, new types can be added to the end
65  * these *MUST* correspond to the arrays in hackblock.cpp
66  */
67 enum
68 {
69     ebtsBONDS,
70     ebtsANGLES,
71     ebtsPDIHS,
72     ebtsIDIHS,
73     ebtsEXCLS,
74     ebtsCMAP,
75     ebtsNR
76 };
77 //! Names for interaction type entries
78 extern const char* btsNames[ebtsNR];
79 //! Numbers for atoms in the interactions.
80 extern const int btsNiatoms[ebtsNR];
81
82 /* if changing any of these structs, make sure that all of the
83    free/clear/copy/merge_t_* functions stay updated */
84
85 /*! \libinternal \brief
86  * Information about single bonded interaction.
87  */
88 struct BondedInteraction
89 {
90     //! Atom names in the bond.
91     std::array<std::string, MAXATOMLIST> a;
92     /*! \brief
93      * Optional define string which gets copied from
94      * .rtp/.tdb to .top and will be parsed by cpp
95      * during grompp.
96      */
97     std::string s;
98     //! Has the entry been found?
99     bool match = false;
100     //! Get name of first atom in bonded interaction.
101     const std::string& ai() const { return a[0]; }
102     //! Get name of second atom in bonded interaction..
103     const std::string& aj() const { return a[1]; }
104     //! Get name of third atom in bonded interaction.
105     const std::string& ak() const { return a[2]; }
106     //! Get name of fourth atom in bonded interaction.
107     const std::string& al() const { return a[3]; }
108     //! Get name of fifth atom in bonded interaction.
109     const std::string& am() const { return a[4]; }
110 };
111
112 /*! \libinternal \brief
113  * Accumulation of different bonded types for preprocessing.
114  * \todo This should be merged with BondedInteraction.
115  */
116 struct BondedInteractionList
117 {
118     //! The type of bonded interaction.
119     int type = -1;
120     //! The actual bonded interactions.
121     std::vector<BondedInteraction> b;
122 };
123
124 /*! \libinternal \brief
125  * Information about preprocessing residues.
126  */
127 struct PreprocessResidue
128 {
129     //! Name of the residue.
130     std::string resname;
131     //! The base file name this rtp entry was read from.
132     std::string filebase;
133     //! Atom data.
134     std::vector<t_atom> atom;
135     //! Atom names.
136     std::vector<char**> atomname;
137     //! Charge group numbers.
138     std::vector<int> cgnr;
139     //! Delete autogenerated dihedrals or not.
140     bool bKeepAllGeneratedDihedrals = false;
141     //! Number of bonded exclusions.
142     int nrexcl = -1;
143     //! If Hydrogen only 1-4 interactions should be generated.
144     bool bGenerateHH14Interactions = false;
145     //! Delete dihedrals also defined by impropers.
146     bool bRemoveDihedralIfWithImproper = false;
147     //! List of bonded interactions to potentially add.
148     std::array<BondedInteractionList, ebtsNR> rb;
149     //! Get number of atoms in residue.
150     int natom() const { return atom.size(); }
151 };
152
153 //! Declare different types of hacks for later check.
154 enum class MoleculePatchType
155 {
156     //! Hack adds atom to structure/rtp.
157     Add,
158     //! Hack deletes atom.
159     Delete,
160     //! Hack replaces atom.
161     Replace
162 };
163
164 /*! \libinternal \brief
165  * Block to modify individual residues
166  */
167 struct MoleculePatch
168 {
169     //! Number of new are deleted atoms. NOT always equal to atom.size()!
170     int nr;
171     //! Old name for entry.
172     std::string oname;
173     //! New name for entry.
174     std::string nname;
175     //! New atom data.
176     std::vector<t_atom> atom;
177     //! Chargegroup number.
178     int cgnr = NOTSET;
179     //! Type of attachement.
180     int tp = 0;
181     //! Number of control atoms.
182     int nctl = 0;
183     //! Name of control atoms.
184     std::array<std::string, 4> a;
185     //! Is an atom to be hacked already present?
186     bool bAlreadyPresent = false;
187     //! Are coordinates for a new atom already set?
188     bool bXSet = false;
189     //! New position for hacked atom.
190     rvec newx = { NOTSET };
191     //! New atom index number after additions.
192     int newi = -1;
193
194     /*! \brief
195      * Get type of hack.
196      *
197      * This depends on the setting of oname and nname
198      * for legacy reasons. If oname is empty, we are adding,
199      * if oname is set and nname is empty, an atom is deleted,
200      * if both are set replacement is going on. If both are unset,
201      * an error is thrown.
202      */
203     MoleculePatchType type() const;
204
205     //! Control atom i name.
206     const std::string& ai() const { return a[0]; }
207     //! Control atom j name.
208     const std::string& aj() const { return a[1]; }
209     //! Control atom k name.
210     const std::string& ak() const { return a[2]; }
211     //! Control atom l name.
212     const std::string& al() const { return a[3]; }
213 };
214 /*! \libinternal \brief
215  * A set of modifications to apply to atoms.
216  */
217 struct MoleculePatchDatabase
218 {
219     //! Name of block
220     std::string name;
221     //! File that entry was read from.
222     std::string filebase;
223     //! List of changes to atoms.
224     std::vector<MoleculePatch> hack;
225     //! List of bonded interactions to potentially add.
226     std::array<BondedInteractionList, ebtsNR> rb;
227     //! Number of atoms to modify
228     int nhack() const { return hack.size(); }
229 };
230
231 /*!\brief
232  * Reset modification block.
233  *
234  * \param[inout] globalPatches Block to reset.
235  * \todo Remove once constructor/destructor takes care of all of this.
236  */
237 void clearModificationBlock(MoleculePatchDatabase* globalPatches);
238
239 /*!\brief
240  * Copy residue information.
241  *
242  * \param[in] s Source information.
243  * \param[in] d Destination to copy to.
244  * \param[inout] symtab Symbol table for names.
245  * \todo Remove once copy can be done directly.
246  */
247 void copyPreprocessResidues(const PreprocessResidue& s, PreprocessResidue* d, t_symtab* symtab);
248
249 /*! \brief
250  * Add bond information in \p s to \p d.
251  *
252  * \param[in] s Source information to copy.
253  * \param[inout] d Destination to copy to.
254  * \param[in] bMin don't copy bondeds with atoms starting with '-'.
255  * \param[in] bPlus don't copy bondeds with atoms starting with '+'.
256  * \returns if bonds were removed at the termini.
257  */
258 bool mergeBondedInteractionList(gmx::ArrayRef<const BondedInteractionList> s,
259                                 gmx::ArrayRef<BondedInteractionList>       d,
260                                 bool                                       bMin,
261                                 bool                                       bPlus);
262
263 /*! \brief
264  * Copy all information from datastructure.
265  *
266  * \param[in] s Source information.
267  * \param[inout] d Destination to copy to.
268  */
269 void copyModificationBlocks(const MoleculePatchDatabase& s, MoleculePatchDatabase* d);
270
271 /*!\brief
272  * Add the individual modifications in \p s to \p d.
273  *
274  * \param[in] s Source information.
275  * \param[inout] d Destination to copy to.
276  */
277 void mergeAtomModifications(const MoleculePatchDatabase& s, MoleculePatchDatabase* d);
278
279 //! \copydoc mergeAtomModifications
280 void mergeAtomAndBondModifications(const MoleculePatchDatabase& s, MoleculePatchDatabase* d);
281
282 #endif