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,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.
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) :
75 atomName(aName), residueName(rName), isAvailable(false), value(0.0)
81 std::string residueName;
82 //! Is property available.
84 //! Value set for property.
88 //! Conglomeration of atom property entries.
91 //! Has property been set.
93 //! Database the property is coming from.
95 //! Default value for property.
97 //! Basic entries for properties.
98 std::vector<BaseEntry> entry;
101 //! Implementation detail type for Atomproperties.
102 class AtomProperties::Impl
105 //! Should user be warned about error.
106 bool bWarned = false;
107 //! Should user be warned about vdW not found.
108 bool bWarnVDW = false;
109 //! The different atom properties.
110 AtomProperty prop[epropNR];
111 //! The residue types.
112 ResidueTypeMap residueTypeMap = residueTypeMapFromLibraryFile("residuetypes.dat");
116 * Find number of matching characters in entry.
118 * If not all characters are matching, return NOTFOUND.
119 * If the length of the database entry is different from the search,
120 * also return NOTFOUND.
122 * \param[in] search Entry to compare to database.
123 * \param[in] database Name of the database entry to compare to.
124 * \returns Number of matching characters or NOTFOUND.
126 static int compareToDatabase(const std::string& search, const std::string& database)
128 if (database.length() > search.length())
133 for (size_t i = 0; i < database.length(); i++)
135 if (search[i] == database[i])
140 if (matches == database.length())
151 * Finds the index for the property being searched.
153 * \param[in] ap Property to search for.
154 * \param[in] residueTypeMap Residuetypes in database.
155 * \param[in] residueName The name of the residue to look for.
156 * \param[in] atomName The name of the atom to look for.
157 * \param[in] bExact Do we have the correct match.
158 * \returns The index for the property.
160 static int findPropertyIndex(AtomProperty* ap,
161 const ResidueTypeMap& residueTypeMap,
162 const std::string& residueName,
163 const std::string& atomName,
168 bool bProtein = namedResidueHasType(residueTypeMap, residueName, "Protein");
169 bool bProtWild = residueName == "AAA";
170 int malen = NOTFOUND;
171 int mrlen = NOTFOUND;
172 for (size_t i = 0; (i < ap->entry.size()); i++)
174 int rlen = compareToDatabase(residueName, ap->entry[i].residueName);
175 if (rlen == NOTFOUND)
177 if ((ap->entry[i].residueName == "*") || (ap->entry[i].residueName == "???"))
181 else if (ap->entry[i].residueName == "AAA")
186 int alen = compareToDatabase(atomName, ap->entry[i].atomName);
187 if ((alen > NOTFOUND) && (rlen > NOTFOUND))
189 if (((alen > malen) && (rlen >= mrlen)) || ((rlen > mrlen) && (alen >= malen)))
198 *bExact = ((malen == static_cast<long int>(atomName.length()))
199 && ((mrlen == static_cast<long int>(residueName.length())) || ((mrlen == WILDPROT) && bProtWild)
200 || ((mrlen == WILDCARD) && !bProtein && !bProtWild)));
204 fprintf(debug, "searching residue: %4s atom: %4s\n", residueName.c_str(), atomName.c_str());
207 fprintf(debug, " not successful\n");
213 ap->entry[j].residueName.c_str(),
214 ap->entry[j].atomName.c_str());
221 * Add new property to list.
223 * \param[in] ap Atomproperty to add.
224 * \param[in] residueTypeMap Residue type database to use.
225 * \param[in] residueName Name of the residue.
226 * \param[in] atomName Name of the atom.
227 * \param[in] propValue Value of property.
228 * \param[in] line Where to add property.
230 static void addProperty(AtomProperty* ap,
231 const ResidueTypeMap& residueTypeMap,
232 const std::string& residueName,
233 const std::string& atomName,
238 int j = findPropertyIndex(ap, residueTypeMap, residueName, atomName, &bExact);
242 ap->entry.emplace_back(BaseEntry(atomName, residueName));
244 j = ap->entry.size() - 1;
246 if (ap->entry[j].isAvailable)
248 if (ap->entry[j].value == propValue)
251 "Warning double identical entries for %s %s %g on line %d in file %s\n",
261 "Warning double different entries %s %s %g and %g on line %d in file %s\n"
262 "Using last entry (%g)\n",
270 ap->entry[j].value = propValue;
275 ap->entry[j].isAvailable = TRUE;
276 ap->entry[j].value = propValue;
281 * Read property value into structure.
283 * \param[in] ap Atomproperty to be read in.
284 * \param[in] residueTypeMap Library of residue types.
285 * \param[in] factor Scaling factor for property.
287 static void readProperty(AtomProperty* ap, const ResidueTypeMap& residueTypeMap, double factor)
289 char line[STRLEN], resnm[32], atomnm[32];
291 gmx::FilePtr fp = gmx::openLibraryFile(ap->db);
293 while (get_a_line(fp.get(), line, STRLEN))
297 if (sscanf(line, "%31s %31s %20lf", resnm, atomnm, &pp) == 3)
300 addProperty(ap, residueTypeMap, resnm, atomnm, pp, line_no);
304 fprintf(stderr, "WARNING: Error in file %s at line %d ignored\n", ap->db.c_str(), line_no);
311 * Set value for properties.
313 * \param[in] ap Atomproperty to set.
314 * \param[in] residueTypeMap Library of residue types.
315 * \param[in] eprop Which property to set.
316 * \param[in] haveBeenWarned If we already set a warning before
317 * \returns True of warning should be printed.
319 static bool setProperties(AtomProperty* ap, const ResidueTypeMap& residueTypeMap, int eprop, bool haveBeenWarned)
321 const char* fns[epropNR] = {
322 "atommass.dat", "vdwradii.dat", "dgsolv.dat", "electroneg.dat", "elements.dat"
324 double fac[epropNR] = { 1.0, 1.0, 418.4, 1.0, 1.0 };
325 double def[epropNR] = { 12.011, 0.14, 0.0, 2.2, -1 };
327 bool printWarning = false;
331 ap->def = def[eprop];
332 readProperty(ap, residueTypeMap, fac[eprop]);
336 fprintf(debug, "Entries in %s: %zu\n", ap->db.c_str(), ap->entry.size());
339 if ((!haveBeenWarned && (eprop == epropMass)) || (eprop == epropVDW))
347 AtomProperties::AtomProperties() : impl_(new Impl) {}
349 AtomProperties::~AtomProperties() {}
351 AtomProperty* AtomProperties::prop(int eprop)
353 return &impl_->prop[eprop];
356 //! Print warning that vdW radii and masses are guessed.
357 static void printWarning()
360 "WARNING: Masses and atomic (Van der Waals) radii will be guessed\n"
361 " based on residue and atom names, since they could not be\n"
362 " definitively assigned from the information in your input\n"
363 " files. These guessed numbers might deviate from the mass\n"
364 " and radius of the atom type. Please check the output\n"
365 " files if necessary.\n\n");
368 static void printvdwWarning(FILE* fp)
373 "NOTE: From version 5.0 %s uses the Van der Waals radii\n",
374 gmx::getProgramContext().displayName());
375 fprintf(fp, "from the source below. This means the results may be different\n");
376 fprintf(fp, "compared to previous GROMACS versions.\n");
377 please_cite(fp, "Bondi1964a");
381 bool AtomProperties::setAtomProperty(int eprop,
382 const std::string& residueName,
383 const std::string& atomName,
386 std::string tmpAtomName, tmpResidueName;
389 if (setProperties(prop(eprop), impl_->residueTypeMap, eprop, impl_->bWarned))
392 impl_->bWarned = true;
394 if (isdigit(atomName[0]))
396 /* put digit after atomname */
397 tmpAtomName.append(atomName.substr(1));
398 tmpAtomName.append(1, atomName[0]);
402 tmpAtomName = atomName;
404 const int j = findPropertyIndex(
405 &(impl_->prop[eprop]), impl_->residueTypeMap, residueName, tmpAtomName, &bExact);
407 if (eprop == epropVDW && !impl_->bWarnVDW)
409 printvdwWarning(stdout);
410 impl_->bWarnVDW = true;
414 *value = impl_->prop[eprop].entry[j].value;
419 *value = impl_->prop[eprop].def;
425 std::string AtomProperties::elementFromAtomNumber(int atomNumber)
427 if (setProperties(prop(epropElement), impl_->residueTypeMap, epropElement, impl_->bWarned))
430 impl_->bWarned = true;
432 for (const auto& e : prop(epropElement)->entry)
434 if (std::round(e.value) == atomNumber)
442 int AtomProperties::atomNumberFromElement(const char* element)
444 if (setProperties(prop(epropElement), impl_->residueTypeMap, epropElement, impl_->bWarned))
447 impl_->bWarned = true;
449 for (const auto& e : prop(epropElement)->entry)
451 if (gmx_strcasecmp(e.atomName.c_str(), element) == 0)
453 return gmx::roundToInt(e.value);