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