SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / topology / residuetypes.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 #include "gmxpre.h"
37
38 #include "residuetypes.h"
39
40 #include <string>
41
42 #include "gromacs/utility/cstringutil.h"
43 #include "gromacs/utility/fatalerror.h"
44 #include "gromacs/utility/futil.h"
45 #include "gromacs/utility/strdb.h"
46
47 //! Definition for residue type that is not known.
48 const ResidueType c_undefinedResidueType = "Other";
49
50 void addResidue(ResidueTypeMap* residueTypeMap, const ResidueName& residueName, const ResidueType& residueType)
51 {
52     if (auto [foundIt, insertionTookPlace] = residueTypeMap->insert({ residueName, residueType });
53         !insertionTookPlace)
54     {
55         if (!gmx::equalCaseInsensitive(foundIt->second, residueType))
56         {
57             fprintf(stderr,
58                     "Warning: Residue '%s' already present with type '%s' in database, ignoring "
59                     "new type '%s'.\n",
60                     residueName.c_str(),
61                     foundIt->second.c_str(),
62                     residueType.c_str());
63         }
64     }
65 }
66
67 ResidueTypeMap residueTypeMapFromLibraryFile(const std::string& residueTypesDatFilename)
68 {
69     char line[STRLEN];
70     char resname[STRLEN], restype[STRLEN], dum[STRLEN];
71
72     gmx::FilePtr db = gmx::openLibraryFile(residueTypesDatFilename);
73
74     ResidueTypeMap residueTypeMap;
75     while (get_a_line(db.get(), line, STRLEN))
76     {
77         strip_comment(line);
78         trim(line);
79         if (line[0] != '\0')
80         {
81             if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
82             {
83                 gmx_fatal(
84                         FARGS,
85                         "Incorrect number of columns (2 expected) for line in residuetypes.dat  ");
86             }
87             addResidue(&residueTypeMap, resname, restype);
88         }
89     }
90     return residueTypeMap;
91 }
92
93 bool namedResidueHasType(const ResidueTypeMap& residueTypeMap,
94                          const ResidueName&    residueName,
95                          const ResidueType&    residueType)
96 {
97     if (auto foundIt = residueTypeMap.find(residueName); foundIt != residueTypeMap.end())
98     {
99         return gmx::equalCaseInsensitive(residueType, foundIt->second);
100     }
101     return false;
102 }
103
104 ResidueType typeOfNamedDatabaseResidue(const ResidueTypeMap& residueTypeMap, const ResidueName& residueName)
105 {
106     if (auto foundIt = residueTypeMap.find(residueName); foundIt != residueTypeMap.end())
107     {
108         return foundIt->second;
109     }
110     return c_undefinedResidueType;
111 }