clang-tidy-11 fixes for gmxpreprocess
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / hackblock.cpp
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,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020,2021, 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 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "hackblock.h"
42
43 #include <cstring>
44
45 #include "gromacs/gmxpreprocess/notset.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/topology/atoms.h"
49 #include "gromacs/topology/symtab.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/enumerationhelpers.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/smalloc.h"
57 #include "gromacs/utility/stringcompare.h"
58
59 const char* enumValueToString(BondedTypes enumValue)
60 {
61     /* these MUST correspond to the enum in hackblock.h */
62     constexpr gmx::EnumerationArray<BondedTypes, const char*> bondedTypeNames = {
63         "bonds", "angles", "dihedrals", "impropers", "exclusions", "cmap"
64     };
65     return bondedTypeNames[enumValue];
66 }
67
68 int enumValueToNumIAtoms(BondedTypes enumValue)
69 {
70     constexpr gmx::EnumerationArray<BondedTypes, int> bondedTypeIAtoms = { 2, 3, 4, 4, 2, 5 };
71     return bondedTypeIAtoms[enumValue];
72 }
73
74 MoleculePatchType MoleculePatch::type() const
75 {
76     if (oname.empty() && !nname.empty())
77     {
78         return MoleculePatchType::Add;
79     }
80     else if (!oname.empty() && nname.empty())
81     {
82         return MoleculePatchType::Delete;
83     }
84     else if (!oname.empty() && !nname.empty())
85     {
86         return MoleculePatchType::Replace;
87     }
88     else
89     {
90         GMX_THROW(gmx::InvalidInputError("Unknown type of atom modification"));
91     }
92 }
93
94 void clearModificationBlock(MoleculePatchDatabase* globalPatches)
95 {
96     globalPatches->name.clear();
97     globalPatches->hack.clear();
98     for (auto bondedsList : globalPatches->rb)
99     {
100         bondedsList.b.clear();
101     }
102 }
103
104 #define safe_strdup(str) (((str) != NULL) ? gmx_strdup(str) : NULL)
105
106 static bool contains_char(const BondedInteraction& s, char c)
107 {
108
109     bool bRet = false;
110     for (int i = 0; i < MAXATOMLIST; i++)
111     {
112         if (!s.a[i].empty() && s.a[i][0] == c)
113         {
114             bRet = true;
115         }
116     }
117
118     return bRet;
119 }
120
121 static int rbonded_find_atoms_in_list(const BondedInteraction&               b,
122                                       gmx::ArrayRef<const BondedInteraction> blist,
123                                       int                                    natoms)
124 {
125     int foundPos = -1;
126
127     for (auto it = blist.begin(); (it != blist.end()) && (foundPos < 0); it++)
128     {
129         bool atomsMatch = true;
130         for (int k = 0; k < natoms && atomsMatch; k++)
131         {
132             atomsMatch = atomsMatch && (b.a[k] == it->a[k]);
133         }
134         /* Try reverse if forward match did not work */
135         if (!atomsMatch)
136         {
137             atomsMatch = true;
138             for (int k = 0; k < natoms && atomsMatch; k++)
139             {
140                 atomsMatch = atomsMatch && (b.a[k] == it->a[natoms - 1 - k]);
141             }
142         }
143         if (atomsMatch)
144         {
145             foundPos = std::distance(blist.begin(), it);
146             /* If all the atoms AND all the parameters match, it is likely that
147              * the user made a copy-and-paste mistake (since it would be much cheaper
148              * to just bump the force constant 2x if you really want it twice).
149              * Since we only have the unparsed string here we can only detect
150              * EXACT matches (including identical whitespace).
151              */
152             if (b.s == it->s)
153             {
154                 gmx_warning("Duplicate line found in or between hackblock and rtp entries");
155             }
156         }
157     }
158     return foundPos;
159 }
160
161 bool mergeBondedInteractionList(gmx::ArrayRef<const BondedInteractionList> s,
162                                 gmx::ArrayRef<BondedInteractionList>       d,
163                                 bool                                       bMin,
164                                 bool                                       bPlus)
165 {
166     bool bBondsRemoved = false;
167     for (auto i : gmx::EnumerationWrapper<BondedTypes>{})
168     {
169         int value = static_cast<int>(i);
170         if (!s[value].b.empty())
171         {
172             /* Record how many bonds we have in the destination when we start.
173              *
174              * If an entry is present in the hackblock (destination), we will
175              * not add the one from the main rtp, since the point is for hackblocks
176              * to overwrite it. However, if there is no hackblock entry we do
177              * allow multiple main rtp entries since some forcefield insist on that.
178              *
179              * We accomplish this by checking the position we find an entry in,
180              * rather than merely checking whether it exists at all.
181              * If that index is larger than the original (hackblock) destination
182              * size, it was added from the main rtp, and then we will allow more
183              * such entries. In contrast, if the entry found has a lower index
184              * it is a hackblock entry meant to override the main rtp, and then
185              * we don't add the main rtp one.
186              */
187             int nbHackblockStart = d[value].b.size();
188
189             for (const auto& b : s[value].b)
190             {
191                 /* Check if this bonded string already exists before adding.
192                  * We are merging from the main RTP to the hackblocks, so this
193                  * will mean the hackblocks overwrite the man RTP, as intended.
194                  */
195                 int index = rbonded_find_atoms_in_list(b, d[value].b, enumValueToNumIAtoms(i));
196                 /* - If we did not find this interaction at all, the index will be -1,
197                  *   and then we should definitely add it to the merged hackblock and rtp.
198                  *
199                  * Alternatively, if it was found, index will be >=0.
200                  * - In case this index is lower than the original number of entries,
201                  *   it is already present as a *hackblock* entry, and those should
202                  *   always override whatever we have listed in the RTP. Thus, we
203                  *   should just keep that one and not add anything from the RTP.
204                  * - Finally, if it was found, but with an index higher than
205                  *   the original number of entries, it comes from the RTP rather
206                  *   than hackblock, and then we must have added it ourselves
207                  *   in a previous iteration. In that case it is a matter of
208                  *   several entries for the same sequence of atoms, and we allow
209                  *   that in the RTP. In this case we should simply copy all of
210                  *   them, including this one.
211                  */
212                 if (index < 0 || index >= nbHackblockStart)
213                 {
214                     if (!(bMin && contains_char(b, '-')) && !(bPlus && contains_char(b, '+')))
215                     {
216                         d[value].b.push_back(b);
217                     }
218                     else if (i == BondedTypes::Bonds)
219                     {
220                         bBondsRemoved = true;
221                     }
222                 }
223                 else
224                 {
225                     /* This is the common case where a hackblock entry simply
226                      * overrides the RTP, so we cannot warn here.
227                      */
228                 }
229             }
230         }
231     }
232     return bBondsRemoved;
233 }
234
235 void copyPreprocessResidues(const PreprocessResidue& s, PreprocessResidue* d, t_symtab* symtab)
236 {
237     *d = s;
238     d->atom.clear();
239     for (const auto& a : s.atom)
240     {
241         d->atom.push_back(a);
242     }
243     d->atomname.clear();
244     for (const auto& a : s.atomname)
245     {
246         d->atomname.push_back(put_symtab(symtab, *a));
247     }
248     d->cgnr.clear();
249     for (const auto& c : s.cgnr)
250     {
251         d->cgnr.push_back(c);
252     }
253     for (auto i : gmx::EnumerationWrapper<BondedTypes>{})
254     {
255         d->rb[i].type = s.rb[i].type;
256         d->rb[i].b.clear();
257     }
258     mergeBondedInteractionList(s.rb, d->rb, FALSE, FALSE);
259 }
260
261 void mergeAtomModifications(const MoleculePatchDatabase& s, MoleculePatchDatabase* d)
262 {
263     for (const auto& h : s.hack)
264     {
265         d->hack.push_back(h);
266     }
267 }
268
269 void mergeAtomAndBondModifications(const MoleculePatchDatabase& s, MoleculePatchDatabase* d)
270 {
271     mergeAtomModifications(s, d);
272     mergeBondedInteractionList(s.rb, d->rb, FALSE, FALSE);
273 }
274
275 void copyModificationBlocks(const MoleculePatchDatabase& s, MoleculePatchDatabase* d)
276 {
277     *d      = s;
278     d->name = s.name;
279     d->hack.clear();
280     for (auto bondedList : d->rb)
281     {
282         bondedList.b.clear();
283     }
284     mergeAtomAndBondModifications(s, d);
285 }
286
287 #undef safe_strdup