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