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