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, 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.
50 #include "gromacs/gmxpreprocess/fflibutil.h"
51 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
52 #include "gromacs/gmxpreprocess/grompp_impl.h"
53 #include "gromacs/gmxpreprocess/notset.h"
54 #include "gromacs/gmxpreprocess/pgutil.h"
55 #include "gromacs/topology/atoms.h"
56 #include "gromacs/topology/symtab.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/strdb.h"
63 #include "hackblock.h"
65 PreprocessingAtomTypes read_atype(const char* ffdir, t_symtab* tab)
68 char buf[STRLEN], name[STRLEN];
72 std::vector<std::string> files = fflib_search_file_end(ffdir, ".atp", TRUE);
74 PreprocessingAtomTypes at;
76 for (const auto& filename : files)
78 in = fflib_open(filename);
81 /* Skip blank or comment-only lines */
84 if (fgets2(buf, STRLEN, in) != nullptr)
89 } while ((feof(in) == 0) && strlen(buf) == 0);
91 if (sscanf(buf, "%s%lf", name, &m) == 2)
94 at.addType(tab, *a, name, InteractionOfType({}, {}), 0, 0);
99 gmx_fatal(FARGS, "Invalid atomtype format: '%s'", buf);
108 static void print_resatoms(FILE* out, const PreprocessingAtomTypes& atype, const PreprocessResidue& rtpDBEntry)
110 fprintf(out, "[ %s ]\n", rtpDBEntry.resname.c_str());
111 fprintf(out, " [ atoms ]\n");
113 for (int j = 0; (j < rtpDBEntry.natom()); j++)
115 int tp = rtpDBEntry.atom[j].type;
116 const char* tpnm = atype.atomNameFromAtomType(tp);
119 gmx_fatal(FARGS, "Incorrect atomtype (%d)", tp);
121 fprintf(out, "%6s %6s %8.3f %6d\n", *(rtpDBEntry.atomname[j]), tpnm,
122 rtpDBEntry.atom[j].q, rtpDBEntry.cgnr[j]);
126 static bool read_atoms(FILE* in, char* line, PreprocessResidue* r0, t_symtab* tab, PreprocessingAtomTypes* atype)
129 char buf[256], buf1[256];
134 r0->atomname.clear();
136 while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
138 if (sscanf(line, "%s%s%lf%d", buf, buf1, &q, &cg) != 4)
142 r0->atomname.push_back(put_symtab(tab, buf));
143 r0->atom.emplace_back();
144 r0->atom.back().q = q;
145 r0->cgnr.push_back(cg);
146 int j = atype->atomTypeFromName(buf1);
150 "Atom type %s (residue %s) not found in atomtype "
152 buf1, r0->resname.c_str());
154 r0->atom.back().type = j;
155 r0->atom.back().m = atype->atomMassFromAtomType(j);
161 static bool read_bondeds(int bt, FILE* in, char* line, PreprocessResidue* rtpDBEntry)
165 while (get_a_line(in, line, STRLEN) && (strchr(line, '[') == nullptr))
169 rtpDBEntry->rb[bt].b.emplace_back();
170 BondedInteraction* newBond = &rtpDBEntry->rb[bt].b.back();
171 for (int j = 0; j < btsNiatoms[bt]; j++)
173 if (sscanf(line + n, "%s%n", str, &ni) == 1)
183 while (isspace(line[n]))
188 newBond->s = line + n;
194 static void print_resbondeds(FILE* out, int bt, const PreprocessResidue& rtpDBEntry)
196 if (!rtpDBEntry.rb[bt].b.empty())
198 fprintf(out, " [ %s ]\n", btsNames[bt]);
200 for (const auto& b : rtpDBEntry.rb[bt].b)
202 for (int j = 0; j < btsNiatoms[bt]; j++)
204 fprintf(out, "%6s ", b.a[j].c_str());
208 fprintf(out, " %s", b.s.c_str());
215 static void check_rtp(gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const std::string& libfn)
217 /* check for double entries, assuming list is already sorted */
218 for (auto it = rtpDBEntry.begin() + 1; it != rtpDBEntry.end(); it++)
221 if (gmx::equalCaseInsensitive(prev->resname, it->resname))
223 fprintf(stderr, "WARNING double entry %s in file %s\n", it->resname.c_str(), libfn.c_str());
228 static int get_bt(char* header)
232 for (i = 0; i < ebtsNR; i++)
234 if (gmx_strcasecmp(btsNames[i], header) == 0)
242 /* print all the ebtsNR type numbers */
243 static void print_resall_header(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
245 fprintf(out, "[ bondedtypes ]\n");
247 "; bonds angles dihedrals impropers all_dihedrals nr_exclusions HH14 "
249 fprintf(out, " %5d %6d %9d %9d %14d %14d %14d %14d\n\n", rtpDBEntry[0].rb[0].type,
250 rtpDBEntry[0].rb[1].type, rtpDBEntry[0].rb[2].type, rtpDBEntry[0].rb[3].type,
251 static_cast<int>(rtpDBEntry[0].bKeepAllGeneratedDihedrals), rtpDBEntry[0].nrexcl,
252 static_cast<int>(rtpDBEntry[0].bGenerateHH14Interactions),
253 static_cast<int>(rtpDBEntry[0].bRemoveDihedralIfWithImproper));
256 void print_resall(FILE* out, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry, const PreprocessingAtomTypes& atype)
258 if (rtpDBEntry.empty())
263 print_resall_header(out, rtpDBEntry);
265 for (const auto& r : rtpDBEntry)
269 print_resatoms(out, atype, r);
270 for (int bt = 0; bt < ebtsNR; bt++)
272 print_resbondeds(out, bt, r);
278 void readResidueDatabase(const std::string& rrdb,
279 std::vector<PreprocessResidue>* rtpDBEntry,
280 PreprocessingAtomTypes* atype,
282 bool bAllowOverrideRTP)
285 char filebase[STRLEN], line[STRLEN], header[STRLEN];
287 int dum1, dum2, dum3;
288 bool bNextResidue, bError;
290 fflib_filename_base(rrdb.c_str(), filebase, STRLEN);
292 in = fflib_open(rrdb);
294 PreprocessResidue header_settings;
296 /* these bonded parameters will overwritten be when *
297 * there is a [ bondedtypes ] entry in the .rtp file */
298 header_settings.rb[ebtsBONDS].type = 1; /* normal bonds */
299 header_settings.rb[ebtsANGLES].type = 1; /* normal angles */
300 header_settings.rb[ebtsPDIHS].type = 1; /* normal dihedrals */
301 header_settings.rb[ebtsIDIHS].type = 2; /* normal impropers */
302 header_settings.rb[ebtsEXCLS].type = 1; /* normal exclusions */
303 header_settings.rb[ebtsCMAP].type = 1; /* normal cmap torsions */
305 header_settings.bKeepAllGeneratedDihedrals = FALSE;
306 header_settings.nrexcl = 3;
307 header_settings.bGenerateHH14Interactions = TRUE;
308 header_settings.bRemoveDihedralIfWithImproper = TRUE;
310 /* Column 5 & 6 aren't really bonded types, but we include
311 * them here to avoid introducing a new section:
312 * Column 5 : This controls the generation of dihedrals from the bonding.
313 * All possible dihedrals are generated automatically. A value of
314 * 1 here means that all these are retained. A value of
315 * 0 here requires generated dihedrals be removed if
316 * * there are any dihedrals on the same central atoms
317 * specified in the residue topology, or
318 * * there are other identical generated dihedrals
319 * sharing the same central atoms, or
320 * * there are other generated dihedrals sharing the
321 * same central bond that have fewer hydrogen atoms
322 * Column 6: Number of bonded neighbors to exclude.
323 * Column 7: Generate 1,4 interactions between two hydrogen atoms
324 * Column 8: Remove proper dihedrals if centered on the same bond
325 * as an improper dihedral
327 get_a_line(in, line, STRLEN);
328 if (!get_header(line, header))
330 gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
332 if (gmx::equalCaseInsensitive("bondedtypes", header, 5))
334 get_a_line(in, line, STRLEN);
335 if ((nparam = sscanf(line, "%d %d %d %d %d %d %d %d", &header_settings.rb[ebtsBONDS].type,
336 &header_settings.rb[ebtsANGLES].type,
337 &header_settings.rb[ebtsPDIHS].type, &header_settings.rb[ebtsIDIHS].type,
338 &dum1, &header_settings.nrexcl, &dum2, &dum3))
341 gmx_fatal(FARGS, "need 4 to 8 parameters in the header of .rtp file %s at line:\n%s\n",
344 header_settings.bKeepAllGeneratedDihedrals = (dum1 != 0);
345 header_settings.bGenerateHH14Interactions = (dum2 != 0);
346 header_settings.bRemoveDihedralIfWithImproper = (dum3 != 0);
347 get_a_line(in, line, STRLEN);
350 fprintf(stderr, "Using default: not generating all possible dihedrals\n");
351 header_settings.bKeepAllGeneratedDihedrals = FALSE;
355 fprintf(stderr, "Using default: excluding 3 bonded neighbors\n");
356 header_settings.nrexcl = 3;
360 fprintf(stderr, "Using default: generating 1,4 H--H interactions\n");
361 header_settings.bGenerateHH14Interactions = TRUE;
366 "Using default: removing proper dihedrals found on the same bond as a proper "
368 header_settings.bRemoveDihedralIfWithImproper = TRUE;
374 "Reading .rtp file without '[ bondedtypes ]' directive,\n"
375 "Will proceed as if the entry was:\n");
376 print_resall_header(stderr, gmx::arrayRefFromArray(&header_settings, 1));
378 /* We don't know the current size of rrtp, but simply realloc immediately */
379 auto oldArrayEnd = rtpDBEntry->end();
382 /* Initialise rtp entry structure */
383 rtpDBEntry->push_back(header_settings);
384 PreprocessResidue* res = &rtpDBEntry->back();
385 if (!get_header(line, header))
387 gmx_fatal(FARGS, "in .rtp file at line:\n%s\n", line);
389 res->resname = header;
390 res->filebase = filebase;
392 get_a_line(in, line, STRLEN);
394 bNextResidue = FALSE;
397 if (!get_header(line, header))
406 /* header is an bonded directive */
407 bError = !read_bondeds(bt, in, line, res);
409 else if (gmx::equalCaseInsensitive("atoms", header, 5))
411 /* header is the atoms directive */
412 bError = !read_atoms(in, line, res, tab, atype);
416 /* else header must be a residue name */
422 gmx_fatal(FARGS, "in .rtp file in residue %s at line:\n%s\n", res->resname.c_str(), line);
424 } while ((feof(in) == 0) && !bNextResidue);
426 if (res->natom() == 0)
428 gmx_fatal(FARGS, "No atoms found in .rtp file in residue %s\n", res->resname.c_str());
431 auto found = std::find_if(rtpDBEntry->begin(), rtpDBEntry->end() - 1,
432 [&res](const PreprocessResidue& entry) {
433 return gmx::equalCaseInsensitive(entry.resname, res->resname);
436 if (found != rtpDBEntry->end() - 1)
438 if (found >= oldArrayEnd)
440 gmx_fatal(FARGS, "Found a second entry for '%s' in '%s'", res->resname.c_str(),
443 if (bAllowOverrideRTP)
446 "Found another rtp entry for '%s' in '%s', ignoring this entry and keeping "
447 "the one from '%s.rtp'\n",
448 res->resname.c_str(), rrdb.c_str(), found->filebase.c_str());
449 /* We should free all the data for this entry.
450 * The current code gives a lot of dangling pointers.
452 rtpDBEntry->erase(rtpDBEntry->end() - 1);
457 "Found rtp entries for '%s' in both '%s' and '%s'. If you want the first "
458 "definition to override the second one, set the -rtpo option of pdb2gmx.",
459 res->resname.c_str(), found->filebase.c_str(), rrdb.c_str());
465 std::sort(rtpDBEntry->begin(), rtpDBEntry->end(), [](const PreprocessResidue& a, const PreprocessResidue& b) {
466 return std::lexicographical_compare(
467 a.resname.begin(), a.resname.end(), b.resname.begin(), b.resname.end(),
468 [](const char& c1, const char& c2) { return std::toupper(c1) < std::toupper(c2); });
471 check_rtp(*rtpDBEntry, rrdb);
474 /************************************************************
478 ***********************************************************/
479 static bool is_sign(char c)
481 return (c == '+' || c == '-');
484 /* Compares if the strins match except for a sign at the end
485 * of (only) one of the two.
487 static int neq_str_sign(const char* a1, const char* a2)
491 l1 = static_cast<int>(strlen(a1));
492 l2 = static_cast<int>(strlen(a2));
493 lm = std::min(l1, l2);
495 if (lm >= 1 && ((l1 == l2 + 1 && is_sign(a1[l1 - 1])) || (l2 == l1 + 1 && is_sign(a2[l2 - 1])))
496 && gmx::equalCaseInsensitive(a1, a2, lm))
506 std::string searchResidueDatabase(const std::string& key, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
508 int nbest, best, besti;
513 /* We want to match at least one character */
515 for (auto it = rtpDBEntry.begin(); it != rtpDBEntry.end(); it++)
517 if (gmx::equalCaseInsensitive(key, it->resname))
519 besti = std::distance(rtpDBEntry.begin(), it);
525 /* Allow a mismatch of at most a sign character (with warning) */
526 int n = neq_str_sign(key.c_str(), it->resname.c_str());
527 if (n >= best && n + 1 >= gmx::index(key.length()) && n + 1 >= gmx::index(it->resname.length()))
533 bestbuf = rtpDBEntry[besti].resname;
537 bestbuf.append(" or ");
538 bestbuf.append(it->resname);
545 besti = std::distance(rtpDBEntry.begin(), it);
553 gmx_fatal(FARGS, "Residue '%s' not found in residue topology database, looks a bit like %s",
554 key.c_str(), bestbuf.c_str());
556 else if (besti == -1)
558 gmx_fatal(FARGS, "Residue '%s' not found in residue topology database", key.c_str());
560 if (!gmx::equalCaseInsensitive(rtpDBEntry[besti].resname, key))
563 "\nWARNING: '%s' not found in residue topology database, "
564 "trying to use '%s'\n\n",
565 key.c_str(), rtpDBEntry[besti].resname.c_str());
568 return rtpDBEntry[besti].resname;
571 gmx::ArrayRef<const PreprocessResidue>::const_iterator
572 getDatabaseEntry(const std::string& rtpname, gmx::ArrayRef<const PreprocessResidue> rtpDBEntry)
574 auto found = std::find_if(rtpDBEntry.begin(), rtpDBEntry.end(),
575 [&rtpname](const PreprocessResidue& entry) {
576 return gmx::equalCaseInsensitive(rtpname, entry.resname);
578 if (found == rtpDBEntry.end())
580 /* This should never happen, since searchResidueDatabase should have been called
581 * before calling getDatabaseEntry.
583 gmx_fatal(FARGS, "Residue type '%s' not found in residue topology database", rtpname.c_str());