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) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,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.
51 #include "gromacs/gmxpreprocess/fflibutil.h"
52 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
53 #include "gromacs/gmxpreprocess/grompp_impl.h"
54 #include "gromacs/gmxpreprocess/notset.h"
55 #include "gromacs/gmxpreprocess/pgutil.h"
56 #include "gromacs/topology/atoms.h"
57 #include "gromacs/topology/symtab.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/logger.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "gromacs/utility/strdb.h"
64 #include "gromacs/utility/stringtoenumvalueconverter.h"
66 #include "hackblock.h"
68 PreprocessingAtomTypes read_atype(const char* ffdir, t_symtab* tab)
71 char buf[STRLEN], name[STRLEN];
75 std::vector<std::string> files = fflib_search_file_end(ffdir, ".atp", TRUE);
77 PreprocessingAtomTypes at;
79 for (const auto& filename : files)
81 in = fflib_open(filename);
84 /* Skip blank or comment-only lines */
87 if (fgets2(buf, STRLEN, in) != nullptr)
92 } while ((feof(in) == 0) && strlen(buf) == 0);
94 if (sscanf(buf, "%s%lf", name, &m) == 2)
97 at.addType(tab, *a, name, InteractionOfType({}, {}), 0, 0);
101 gmx_fatal(FARGS, "Invalid atomtype format: '%s'", buf);
110 static void print_resatoms(FILE* out, const PreprocessingAtomTypes& atype, const PreprocessResidue& rtpDBEntry)
112 fprintf(out, "[ %s ]\n", rtpDBEntry.resname.c_str());
113 fprintf(out, " [ atoms ]\n");
115 for (int j = 0; (j < rtpDBEntry.natom()); j++)
117 int tp = rtpDBEntry.atom[j].type;
118 auto tpnm = atype.atomNameFromAtomType(tp);
119 if (!tpnm.has_value())
121 gmx_fatal(FARGS, "Incorrect atomtype (%d)", tp);
124 "%6s %6s %8.3f %6d\n",
125 *(rtpDBEntry.atomname[j]),
127 rtpDBEntry.atom[j].q,
132 static bool read_atoms(FILE* in, char* line, PreprocessResidue* r0, t_symtab* tab, PreprocessingAtomTypes* atype)
135 char buf[256], buf1[256];
140 r0->atomname.clear();
142 while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
144 if (sscanf(line, "%s%s%lf%d", buf, buf1, &q, &cg) != 4)
148 r0->atomname.push_back(put_symtab(tab, buf));
149 r0->atom.emplace_back();
150 r0->atom.back().q = q;
151 r0->cgnr.push_back(cg);
152 auto j = atype->atomTypeFromName(buf1);
156 "Atom type %s (residue %s) not found in atomtype "
159 r0->resname.c_str());
161 r0->atom.back().type = *j;
162 r0->atom.back().m = *atype->atomMassFromAtomType(*j);
168 static bool read_bondeds(BondedTypes bt, FILE* in, char* line, PreprocessResidue* rtpDBEntry)
172 while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
176 rtpDBEntry->rb[bt].b.emplace_back();
177 BondedInteraction* newBond = &rtpDBEntry->rb[bt].b.back();
178 for (int j = 0; j < enumValueToNumIAtoms(bt); j++)
180 if (sscanf(line + n, "%s%n", str, &ni) == 1)
190 while (isspace(line[n]))
195 newBond->s = line + n;
201 static void print_resbondeds(FILE* out, BondedTypes bt, const PreprocessResidue& rtpDBEntry)
203 if (!rtpDBEntry.rb[bt].b.empty())
205 fprintf(out, " [ %s ]\n", enumValueToString(bt));
207 for (const auto& b : rtpDBEntry.rb[bt].b)
209 for (int j = 0; j < enumValueToNumIAtoms(bt); j++)
211 fprintf(out, "%6s ", b.a[j].c_str());
215 fprintf(out, " %s", b.s.c_str());
222 static void check_rtp(gmx::ArrayRef<const PreprocessResidue> rtpDBEntry,
223 const std::string& libfn,
224 const gmx::MDLogger& logger)
226 /* check for double entries, assuming list is already sorted */
227 for (auto it = rtpDBEntry.begin() + 1; it != rtpDBEntry.end(); it++)
230 if (gmx::equalCaseInsensitive(prev->resname, it->resname))
232 GMX_LOG(logger.warning)
234 .appendTextFormatted("Double entry %s in file %s", it->resname.c_str(), libfn.c_str());
239 static std::optional<BondedTypes> get_bt(char* header)
241 gmx::StringToEnumValueConverter<BondedTypes, enumValueToString> converter;
242 return converter.valueFrom(header);
245 /* print all the BondedTypes type numbers */
246 static void print_resall_header(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
248 fprintf(out, "[ bondedtypes ]\n");
250 "; bonds angles dihedrals impropers all_dihedrals nr_exclusions HH14 "
253 " %5d %6d %9d %9d %14d %14d %14d %14d\n\n",
254 rtpDBEntry[0].rb[0].type,
255 rtpDBEntry[0].rb[1].type,
256 rtpDBEntry[0].rb[2].type,
257 rtpDBEntry[0].rb[3].type,
258 static_cast<int>(rtpDBEntry[0].bKeepAllGeneratedDihedrals),
259 rtpDBEntry[0].nrexcl,
260 static_cast<int>(rtpDBEntry[0].bGenerateHH14Interactions),
261 static_cast<int>(rtpDBEntry[0].bRemoveDihedralIfWithImproper));
265 static void print_resall_log(const gmx::MDLogger& logger, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
267 GMX_LOG(logger.info).asParagraph().appendTextFormatted("[ bondedtypes ]");
270 .appendTextFormatted(
271 "; bonds angles dihedrals impropers all_dihedrals nr_exclusions HH14 "
275 .appendTextFormatted(" %5d %6d %9d %9d %14d %14d %14d %14d",
276 rtpDBEntry[0].rb[0].type,
277 rtpDBEntry[0].rb[1].type,
278 rtpDBEntry[0].rb[2].type,
279 rtpDBEntry[0].rb[3].type,
280 static_cast<int>(rtpDBEntry[0].bKeepAllGeneratedDihedrals),
281 rtpDBEntry[0].nrexcl,
282 static_cast<int>(rtpDBEntry[0].bGenerateHH14Interactions),
283 static_cast<int>(rtpDBEntry[0].bRemoveDihedralIfWithImproper));
287 void print_resall(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const PreprocessingAtomTypes& atype)
289 if (rtpDBEntry.empty())
294 print_resall_header(out, rtpDBEntry);
296 for (const auto& r : rtpDBEntry)
300 print_resatoms(out, atype, r);
301 for (auto bt : gmx::EnumerationWrapper<BondedTypes>{})
303 print_resbondeds(out, bt, r);
309 void readResidueDatabase(const std::string& rrdb,
310 std::vector<PreprocessResidue>* rtpDBEntry,
311 PreprocessingAtomTypes* atype,
313 const gmx::MDLogger& logger,
314 bool bAllowOverrideRTP)
317 char filebase[STRLEN], line[STRLEN], header[STRLEN];
319 int dum1, dum2, dum3;
320 bool bNextResidue, bError;
322 fflib_filename_base(rrdb.c_str(), filebase, STRLEN);
324 in = fflib_open(rrdb);
326 PreprocessResidue header_settings;
328 /* these bonded parameters will overwritten be when *
329 * there is a [ bondedtypes ] entry in the .rtp file */
330 header_settings.rb[BondedTypes::Bonds].type = 1; /* normal bonds */
331 header_settings.rb[BondedTypes::Angles].type = 1; /* normal angles */
332 header_settings.rb[BondedTypes::ProperDihedrals].type = 1; /* normal dihedrals */
333 header_settings.rb[BondedTypes::ImproperDihedrals].type = 2; /* normal impropers */
334 header_settings.rb[BondedTypes::Exclusions].type = 1; /* normal exclusions */
335 header_settings.rb[BondedTypes::Cmap].type = 1; /* normal cmap torsions */
337 header_settings.bKeepAllGeneratedDihedrals = FALSE;
338 header_settings.nrexcl = 3;
339 header_settings.bGenerateHH14Interactions = TRUE;
340 header_settings.bRemoveDihedralIfWithImproper = TRUE;
342 /* Column 5 & 6 aren't really bonded types, but we include
343 * them here to avoid introducing a new section:
344 * Column 5 : This controls the generation of dihedrals from the bonding.
345 * All possible dihedrals are generated automatically. A value of
346 * 1 here means that all these are retained. A value of
347 * 0 here requires generated dihedrals be removed if
348 * * there are any dihedrals on the same central atoms
349 * specified in the residue topology, or
350 * * there are other identical generated dihedrals
351 * sharing the same central atoms, or
352 * * there are other generated dihedrals sharing the
353 * same central bond that have fewer hydrogen atoms
354 * Column 6: Number of bonded neighbors to exclude.
355 * Column 7: Generate 1,4 interactions between two hydrogen atoms
356 * Column 8: Remove proper dihedrals if centered on the same bond
357 * as an improper dihedral
359 get_a_line(in, line, STRLEN);
360 if (!get_header(line, header))
362 gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
364 if (gmx::equalCaseInsensitive("bondedtypes", header, 5))
366 get_a_line(in, line, STRLEN);
367 if ((nparam = sscanf(line,
368 "%d %d %d %d %d %d %d %d",
369 &header_settings.rb[BondedTypes::Bonds].type,
370 &header_settings.rb[BondedTypes::Angles].type,
371 &header_settings.rb[BondedTypes::ProperDihedrals].type,
372 &header_settings.rb[BondedTypes::ImproperDihedrals].type,
374 &header_settings.nrexcl,
380 "need 4 to 8 parameters in the header of .rtp file %s at line:\n%s\n",
384 header_settings.bKeepAllGeneratedDihedrals = (dum1 != 0);
385 header_settings.bGenerateHH14Interactions = (dum2 != 0);
386 header_settings.bRemoveDihedralIfWithImproper = (dum3 != 0);
387 get_a_line(in, line, STRLEN);
392 .appendTextFormatted("Using default: not generating all possible dihedrals");
393 header_settings.bKeepAllGeneratedDihedrals = FALSE;
399 .appendTextFormatted("Using default: excluding 3 bonded neighbors");
400 header_settings.nrexcl = 3;
406 .appendTextFormatted("Using default: generating 1,4 H--H interactions");
407 header_settings.bGenerateHH14Interactions = TRUE;
411 GMX_LOG(logger.warning)
413 .appendTextFormatted(
414 "Using default: removing proper dihedrals found on the same bond as a "
416 header_settings.bRemoveDihedralIfWithImproper = TRUE;
421 GMX_LOG(logger.warning)
423 .appendTextFormatted(
424 "Reading .rtp file without '[ bondedtypes ]' directive, "
425 "Will proceed as if the entry was:");
426 print_resall_log(logger, gmx::arrayRefFromArray(&header_settings, 1));
428 /* We don't know the current size of rrtp, but simply realloc immediately */
429 auto oldArrayEnd = rtpDBEntry->end();
432 /* Initialise rtp entry structure */
433 rtpDBEntry->push_back(header_settings);
434 PreprocessResidue* res = &rtpDBEntry->back();
435 if (!get_header(line, header))
437 gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
439 res->resname = header;
440 res->filebase = filebase;
442 get_a_line(in, line, STRLEN);
444 bNextResidue = FALSE;
447 if (!get_header(line, header))
453 auto bt = get_bt(header);
456 /* header is an bonded directive */
457 bError = !read_bondeds(*bt, in, line, res);
459 else if (gmx::equalCaseInsensitive("atoms", header, 5))
461 /* header is the atoms directive */
462 bError = !read_atoms(in, line, res, tab, atype);
466 /* else header must be a residue name */
472 gmx_fatal(FARGS, "in .rtp file in residue %s at line:\n%s\n", res->resname.c_str(), line);
474 } while ((feof(in) == 0) && !bNextResidue);
476 if (res->natom() == 0)
478 gmx_fatal(FARGS, "No atoms found in .rtp file in residue %s\n", res->resname.c_str());
481 auto found = std::find_if(
482 rtpDBEntry->begin(), rtpDBEntry->end() - 1, [&res](const PreprocessResidue& entry) {
483 return gmx::equalCaseInsensitive(entry.resname, res->resname);
486 if (found != rtpDBEntry->end() - 1)
488 if (found >= oldArrayEnd)
490 gmx_fatal(FARGS, "Found a second entry for '%s' in '%s'", res->resname.c_str(), rrdb.c_str());
492 if (bAllowOverrideRTP)
494 GMX_LOG(logger.warning)
496 .appendTextFormatted(
497 "Found another rtp entry for '%s' in '%s',"
498 " ignoring this entry and keeping the one from '%s.rtp'",
499 res->resname.c_str(),
501 found->filebase.c_str());
502 /* We should free all the data for this entry.
503 * The current code gives a lot of dangling pointers.
505 rtpDBEntry->erase(rtpDBEntry->end() - 1);
510 "Found rtp entries for '%s' in both '%s' and '%s'. If you want the first "
511 "definition to override the second one, set the -rtpo option of pdb2gmx.",
512 res->resname.c_str(),
513 found->filebase.c_str(),
520 std::sort(rtpDBEntry->begin(), rtpDBEntry->end(), [](const PreprocessResidue& a, const PreprocessResidue& b) {
521 return std::lexicographical_compare(
526 [](const char& c1, const char& c2) { return std::toupper(c1) < std::toupper(c2); });
529 check_rtp(*rtpDBEntry, rrdb, logger);
532 /************************************************************
536 ***********************************************************/
537 static bool is_sign(char c)
539 return (c == '+' || c == '-');
542 /* Compares if the strins match except for a sign at the end
543 * of (only) one of the two.
545 static int neq_str_sign(const char* a1, const char* a2)
549 l1 = static_cast<int>(strlen(a1));
550 l2 = static_cast<int>(strlen(a2));
551 lm = std::min(l1, l2);
553 if (lm >= 1 && ((l1 == l2 + 1 && is_sign(a1[l1 - 1])) || (l2 == l1 + 1 && is_sign(a2[l2 - 1])))
554 && gmx::equalCaseInsensitive(a1, a2, lm))
564 std::string searchResidueDatabase(const std::string& key,
565 gmx::ArrayRef<const PreprocessResidue> rtpDBEntry,
566 const gmx::MDLogger& logger)
568 int nbest, best, besti;
573 /* We want to match at least one character */
575 for (auto it = rtpDBEntry.begin(); it != rtpDBEntry.end(); it++)
577 if (gmx::equalCaseInsensitive(key, it->resname))
579 besti = std::distance(rtpDBEntry.begin(), it);
585 /* Allow a mismatch of at most a sign character (with warning) */
586 int n = neq_str_sign(key.c_str(), it->resname.c_str());
587 if (n >= best && n + 1 >= gmx::index(key.length()) && n + 1 >= gmx::index(it->resname.length()))
593 bestbuf = rtpDBEntry[besti].resname;
597 bestbuf.append(" or ");
598 bestbuf.append(it->resname);
605 besti = std::distance(rtpDBEntry.begin(), it);
614 "Residue '%s' not found in residue topology database, looks a bit like %s",
618 else if (besti == -1)
620 gmx_fatal(FARGS, "Residue '%s' not found in residue topology database", key.c_str());
622 if (!gmx::equalCaseInsensitive(rtpDBEntry[besti].resname, key))
624 GMX_LOG(logger.warning)
626 .appendTextFormatted(
627 "'%s' not found in residue topology database, "
628 "trying to use '%s'",
630 rtpDBEntry[besti].resname.c_str());
633 return rtpDBEntry[besti].resname;
636 gmx::ArrayRef<const PreprocessResidue>::const_iterator
637 getDatabaseEntry(const std::string& rtpname, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
639 auto found = std::find_if(
640 rtpDBEntry.begin(), rtpDBEntry.end(), [&rtpname](const PreprocessResidue& entry) {
641 return gmx::equalCaseInsensitive(rtpname, entry.resname);
643 if (found == rtpDBEntry.end())
645 /* This should never happen, since searchResidueDatabase should have been called
646 * before calling getDatabaseEntry.
648 gmx_fatal(FARGS, "Residue type '%s' not found in residue topology database", rtpname.c_str());