2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 /* This file is completely threadsafe - keep it that way! */
41 #include "hackblock.h"
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"
59 const char* enumValueToString(BondedTypes enumValue)
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"
65 return bondedTypeNames[enumValue];
68 int enumValueToNumIAtoms(BondedTypes enumValue)
70 constexpr gmx::EnumerationArray<BondedTypes, int> bondedTypeIAtoms = { 2, 3, 4, 4, 2, 5 };
71 return bondedTypeIAtoms[enumValue];
74 MoleculePatchType MoleculePatch::type() const
76 if (oname.empty() && !nname.empty())
78 return MoleculePatchType::Add;
80 else if (!oname.empty() && nname.empty())
82 return MoleculePatchType::Delete;
84 else if (!oname.empty() && !nname.empty())
86 return MoleculePatchType::Replace;
90 GMX_THROW(gmx::InvalidInputError("Unknown type of atom modification"));
94 void clearModificationBlock(MoleculePatchDatabase* globalPatches)
96 globalPatches->name.clear();
97 globalPatches->hack.clear();
98 for (auto bondedsList : globalPatches->rb)
100 bondedsList.b.clear();
104 static bool contains_char(const BondedInteraction& s, char c)
108 for (int i = 0; i < MAXATOMLIST; i++)
110 if (!s.a[i].empty() && s.a[i][0] == c)
119 static int rbonded_find_atoms_in_list(const BondedInteraction& b,
120 gmx::ArrayRef<const BondedInteraction> blist,
125 for (auto it = blist.begin(); (it != blist.end()) && (foundPos < 0); it++)
127 bool atomsMatch = true;
128 for (int k = 0; k < natoms && atomsMatch; k++)
130 atomsMatch = atomsMatch && (b.a[k] == it->a[k]);
132 /* Try reverse if forward match did not work */
136 for (int k = 0; k < natoms && atomsMatch; k++)
138 atomsMatch = atomsMatch && (b.a[k] == it->a[natoms - 1 - k]);
143 foundPos = std::distance(blist.begin(), it);
144 /* If all the atoms AND all the parameters match, it is likely that
145 * the user made a copy-and-paste mistake (since it would be much cheaper
146 * to just bump the force constant 2x if you really want it twice).
147 * Since we only have the unparsed string here we can only detect
148 * EXACT matches (including identical whitespace).
152 gmx_warning("Duplicate line found in or between hackblock and rtp entries");
159 bool mergeBondedInteractionList(gmx::ArrayRef<const BondedInteractionList> s,
160 gmx::ArrayRef<BondedInteractionList> d,
164 bool bBondsRemoved = false;
165 for (auto i : gmx::EnumerationWrapper<BondedTypes>{})
167 int value = static_cast<int>(i);
168 if (!s[value].b.empty())
170 /* Record how many bonds we have in the destination when we start.
172 * If an entry is present in the hackblock (destination), we will
173 * not add the one from the main rtp, since the point is for hackblocks
174 * to overwrite it. However, if there is no hackblock entry we do
175 * allow multiple main rtp entries since some forcefield insist on that.
177 * We accomplish this by checking the position we find an entry in,
178 * rather than merely checking whether it exists at all.
179 * If that index is larger than the original (hackblock) destination
180 * size, it was added from the main rtp, and then we will allow more
181 * such entries. In contrast, if the entry found has a lower index
182 * it is a hackblock entry meant to override the main rtp, and then
183 * we don't add the main rtp one.
185 int nbHackblockStart = d[value].b.size();
187 for (const auto& b : s[value].b)
189 /* Check if this bonded string already exists before adding.
190 * We are merging from the main RTP to the hackblocks, so this
191 * will mean the hackblocks overwrite the man RTP, as intended.
193 int index = rbonded_find_atoms_in_list(b, d[value].b, enumValueToNumIAtoms(i));
194 /* - If we did not find this interaction at all, the index will be -1,
195 * and then we should definitely add it to the merged hackblock and rtp.
197 * Alternatively, if it was found, index will be >=0.
198 * - In case this index is lower than the original number of entries,
199 * it is already present as a *hackblock* entry, and those should
200 * always override whatever we have listed in the RTP. Thus, we
201 * should just keep that one and not add anything from the RTP.
202 * - Finally, if it was found, but with an index higher than
203 * the original number of entries, it comes from the RTP rather
204 * than hackblock, and then we must have added it ourselves
205 * in a previous iteration. In that case it is a matter of
206 * several entries for the same sequence of atoms, and we allow
207 * that in the RTP. In this case we should simply copy all of
208 * them, including this one.
210 if (index < 0 || index >= nbHackblockStart)
212 if (!(bMin && contains_char(b, '-')) && !(bPlus && contains_char(b, '+')))
214 d[value].b.push_back(b);
216 else if (i == BondedTypes::Bonds)
218 bBondsRemoved = true;
223 /* This is the common case where a hackblock entry simply
224 * overrides the RTP, so we cannot warn here.
230 return bBondsRemoved;
233 void copyPreprocessResidues(const PreprocessResidue& s, PreprocessResidue* d, t_symtab* symtab)
237 for (const auto& a : s.atom)
239 d->atom.push_back(a);
242 for (const auto& a : s.atomname)
244 d->atomname.push_back(put_symtab(symtab, *a));
247 for (const auto& c : s.cgnr)
249 d->cgnr.push_back(c);
251 for (auto i : gmx::EnumerationWrapper<BondedTypes>{})
253 d->rb[i].type = s.rb[i].type;
256 mergeBondedInteractionList(s.rb, d->rb, FALSE, FALSE);
259 void mergeAtomModifications(const MoleculePatchDatabase& s, MoleculePatchDatabase* d)
261 for (const auto& h : s.hack)
263 d->hack.push_back(h);
267 void mergeAtomAndBondModifications(const MoleculePatchDatabase& s, MoleculePatchDatabase* d)
269 mergeAtomModifications(s, d);
270 mergeBondedInteractionList(s.rb, d->rb, FALSE, FALSE);
273 void copyModificationBlocks(const MoleculePatchDatabase& s, MoleculePatchDatabase* d)
278 for (auto bondedList : d->rb)
280 bondedList.b.clear();
282 mergeAtomAndBondModifications(s, d);