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