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,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.
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/math/functions.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/topology/residuetypes.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/pleasecite.h"
57 #include "gromacs/utility/programcontext.h"
58 #include "gromacs/utility/smalloc.h"
59 #include "gromacs/utility/strdb.h"
61 /* NOTFOUND should be smallest, others larger in increasing priority */
70 //! Basic entries in AtomProperty.
73 //! Default constructor.
74 BaseEntry(const std::string& aName, const std::string& rName) :
84 std::string residueName;
85 //! Is property available.
87 //! Value set for property.
91 //! Conglomeration of atom property entries.
94 //! Has property been set.
96 //! Database the property is coming from.
98 //! Default value for property.
100 //! Basic entries for properties.
101 std::vector<BaseEntry> entry;
104 //! Implementation detail type for Atomproperties.
105 class AtomProperties::Impl
108 //! Should user be warned about error.
109 bool bWarned = false;
110 //! Should user be warned about vdW not found.
111 bool bWarnVDW = false;
112 //! The different atom properties.
113 AtomProperty prop[epropNR];
114 //! The residue types.
119 * Find number of matching characters in entry.
121 * If not all characters are matching, return NOTFOUND.
122 * If the length of the database entry is different from the search,
123 * also return NOTFOUND.
125 * \param[in] search Entry to compare to database.
126 * \param[in] database Name of the database entry to compare to.
127 * \returns Number of matching characters or NOTFOUND.
129 static int compareToDatabase(const std::string& search, const std::string& database)
131 if (database.length() > search.length())
136 for (size_t i = 0; i < database.length(); i++)
138 if (search[i] == database[i])
143 if (matches == database.length())
154 * Finds the index for the property being searched.
156 * \param[in] ap Property to search for.
157 * \param[in] restype Residuetypes in database.
158 * \param[in] residueName The name of the residue to look for.
159 * \param[in] atomName The name of the atom to look for.
160 * \param[in] bExact Do we have the correct match.
161 * \returns The index for the property.
163 static int findPropertyIndex(AtomProperty* ap,
164 ResidueType* restype,
165 const std::string& residueName,
166 const std::string& atomName,
171 bool bProtein = restype->namedResidueHasType(residueName, "Protein");
172 bool bProtWild = residueName == "AAA";
173 int malen = NOTFOUND;
174 int mrlen = NOTFOUND;
175 for (size_t i = 0; (i < ap->entry.size()); i++)
177 int rlen = compareToDatabase(residueName, ap->entry[i].residueName);
178 if (rlen == NOTFOUND)
180 if ((ap->entry[i].residueName == "*") || (ap->entry[i].residueName == "???"))
184 else if (ap->entry[i].residueName == "AAA")
189 int alen = compareToDatabase(atomName, ap->entry[i].atomName);
190 if ((alen > NOTFOUND) && (rlen > NOTFOUND))
192 if (((alen > malen) && (rlen >= mrlen)) || ((rlen > mrlen) && (alen >= malen)))
201 *bExact = ((malen == static_cast<long int>(atomName.length()))
202 && ((mrlen == static_cast<long int>(residueName.length())) || ((mrlen == WILDPROT) && bProtWild)
203 || ((mrlen == WILDCARD) && !bProtein && !bProtWild)));
207 fprintf(debug, "searching residue: %4s atom: %4s\n", residueName.c_str(), atomName.c_str());
210 fprintf(debug, " not successful\n");
214 fprintf(debug, " match: %4s %4s\n", ap->entry[j].residueName.c_str(),
215 ap->entry[j].atomName.c_str());
222 * Add new property to list.
224 * \param[in] ap Atomproperty to add.
225 * \param[in] restype Residue type database to use.
226 * \param[in] residueName Name of the residue.
227 * \param[in] atomName Name of the atom.
228 * \param[in] propValue Value of property.
229 * \param[in] line Where to add property.
231 static void addProperty(AtomProperty* ap,
232 ResidueType* restype,
233 const std::string& residueName,
234 const std::string& atomName,
239 int j = findPropertyIndex(ap, restype, residueName, atomName, &bExact);
243 ap->entry.emplace_back(BaseEntry(atomName, residueName));
245 j = ap->entry.size() - 1;
247 if (ap->entry[j].isAvailable)
249 if (ap->entry[j].value == propValue)
251 fprintf(stderr, "Warning double identical entries for %s %s %g on line %d in file %s\n",
252 residueName.c_str(), atomName.c_str(), propValue, line, ap->db.c_str());
257 "Warning double different entries %s %s %g and %g on line %d in file %s\n"
258 "Using last entry (%g)\n",
259 residueName.c_str(), atomName.c_str(), propValue, ap->entry[j].value, line,
260 ap->db.c_str(), propValue);
261 ap->entry[j].value = propValue;
266 ap->entry[j].isAvailable = TRUE;
267 ap->entry[j].value = propValue;
272 * Read property value into structure.
274 * \param[in] ap Atomproperty to be read in.
275 * \param[in] restype Library of residue types.
276 * \param[in] factor Scaling factor for property.
278 static void readProperty(AtomProperty* ap, ResidueType* restype, double factor)
280 char line[STRLEN], resnm[32], atomnm[32];
282 gmx::FilePtr fp = gmx::openLibraryFile(ap->db);
284 while (get_a_line(fp.get(), line, STRLEN))
288 if (sscanf(line, "%31s %31s %20lf", resnm, atomnm, &pp) == 3)
291 addProperty(ap, restype, resnm, atomnm, pp, line_no);
295 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n", ap->db.c_str(), line_no);
302 * Set value for properties.
304 * \param[in] ap Atomproperty to set.
305 * \param[in] restype Library of residue types.
306 * \param[in] eprop Which property to set.
307 * \param[in] haveBeenWarned If we already set a warning before
308 * \returns True of warning should be printed.
310 static bool setProperties(AtomProperty* ap, ResidueType* restype, int eprop, bool haveBeenWarned)
312 const char* fns[epropNR] = { "atommass.dat", "vdwradii.dat", "dgsolv.dat", "electroneg.dat",
314 double fac[epropNR] = { 1.0, 1.0, 418.4, 1.0, 1.0 };
315 double def[epropNR] = { 12.011, 0.14, 0.0, 2.2, -1 };
317 bool printWarning = false;
321 ap->def = def[eprop];
322 readProperty(ap, restype, fac[eprop]);
326 fprintf(debug, "Entries in %s: %zu\n", ap->db.c_str(), ap->entry.size());
329 if ((!haveBeenWarned && (eprop == epropMass)) || (eprop == epropVDW))
337 AtomProperties::AtomProperties() : impl_(new Impl) {}
339 AtomProperties::~AtomProperties() {}
341 AtomProperty* AtomProperties::prop(int eprop)
343 return &impl_->prop[eprop];
346 ResidueType* AtomProperties::restype()
348 return &impl_->restype;
351 //! Print warning that vdW radii and masses are guessed.
352 static void printWarning()
355 "WARNING: Masses and atomic (Van der Waals) radii will be guessed\n"
356 " based on residue and atom names, since they could not be\n"
357 " definitively assigned from the information in your input\n"
358 " files. These guessed numbers might deviate from the mass\n"
359 " and radius of the atom type. Please check the output\n"
360 " files if necessary.\n\n");
363 static void printvdwWarning(FILE* fp)
367 fprintf(fp, "NOTE: From version 5.0 %s uses the Van der Waals radii\n",
368 gmx::getProgramContext().displayName());
369 fprintf(fp, "from the source below. This means the results may be different\n");
370 fprintf(fp, "compared to previous GROMACS versions.\n");
371 please_cite(fp, "Bondi1964a");
375 bool AtomProperties::setAtomProperty(int eprop,
376 const std::string& residueName,
377 const std::string& atomName,
381 std::string tmpAtomName, tmpResidueName;
384 if (setProperties(prop(eprop), restype(), eprop, impl_->bWarned))
387 impl_->bWarned = true;
389 if (isdigit(atomName[0]))
391 /* put digit after atomname */
392 tmpAtomName.append(atomName.substr(1));
393 tmpAtomName.append(1, atomName[0]);
397 tmpAtomName = atomName;
399 j = findPropertyIndex(&(impl_->prop[eprop]), &impl_->restype, residueName, tmpAtomName, &bExact);
401 if (eprop == epropVDW && !impl_->bWarnVDW)
403 printvdwWarning(stdout);
404 impl_->bWarnVDW = true;
408 *value = impl_->prop[eprop].entry[j].value;
413 *value = impl_->prop[eprop].def;
419 std::string AtomProperties::elementFromAtomNumber(int atomNumber)
421 if (setProperties(prop(epropElement), restype(), epropElement, impl_->bWarned))
424 impl_->bWarned = true;
426 for (const auto& e : prop(epropElement)->entry)
428 if (std::round(e.value) == atomNumber)
436 int AtomProperties::atomNumberFromElement(const char* element)
438 if (setProperties(prop(epropElement), restype(), epropElement, impl_->bWarned))
441 impl_->bWarned = true;
443 for (const auto& e : prop(epropElement)->entry)
445 if (gmx_strcasecmp(e.atomName.c_str(), element) == 0)
447 return gmx::roundToInt(e.value);