1584f9f0bec885fe2756676d7e2ef6af3f7e21e2
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / pgutil.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2014,2015,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38
39 #include "gmxpre.h"
40
41 #include "pgutil.h"
42
43 #include <cstring>
44
45 #include "gromacs/topology/atoms.h"
46 #include "gromacs/utility/cstringutil.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/snprintf.h"
49
50 #define BUFSIZE 1024
51 static void atom_not_found(int         fatal_errno,
52                            const char* file,
53                            int         line,
54                            const char* atomname,
55                            int         resind,
56                            const char* resname,
57                            const char* bondtype,
58                            bool        bAllowMissing)
59 {
60     char message_buffer[BUFSIZE];
61     if (strcmp(bondtype, "check") != 0)
62     {
63         if (0 != strcmp(bondtype, "atom"))
64         {
65             snprintf(message_buffer, 1024,
66                      "Residue %d named %s of a molecule in the input file was mapped\n"
67                      "to an entry in the topology database, but the atom %s used in\n"
68                      "an interaction of type %s in that entry is not found in the\n"
69                      "input file. Perhaps your atom and/or residue naming needs to be\n"
70                      "fixed.\n",
71                      resind + 1, resname, atomname, bondtype);
72         }
73         else
74         {
75             snprintf(message_buffer, 1024,
76                      "Residue %d named %s of a molecule in the input file was mapped\n"
77                      "to an entry in the topology database, but the atom %s used in\n"
78                      "that entry is not found in the input file. Perhaps your atom\n"
79                      "and/or residue naming needs to be fixed.\n",
80                      resind + 1, resname, atomname);
81         }
82         if (bAllowMissing)
83         {
84             gmx_warning("WARNING: %s", message_buffer);
85         }
86         else
87         {
88             gmx_fatal(fatal_errno, file, line, "%s", message_buffer);
89         }
90     }
91 }
92
93 int search_atom(const char* type, int start, const t_atoms* atoms, const char* bondtype, bool bAllowMissing)
94 {
95     int           i, resind = -1;
96     bool          bPrevious, bNext;
97     int           natoms = atoms->nr;
98     t_atom*       at     = atoms->atom;
99     char** const* anm    = atoms->atomname;
100
101     bPrevious = (strchr(type, '-') != nullptr);
102     bNext     = (strchr(type, '+') != nullptr);
103
104     if (!bPrevious)
105     {
106         resind = at[start].resind;
107         if (bNext)
108         {
109             /* The next residue */
110             type++;
111             while ((start < natoms) && (at[start].resind == resind))
112             {
113                 start++;
114             }
115             if (start < natoms)
116             {
117                 resind = at[start].resind;
118             }
119         }
120
121         for (i = start; (i < natoms) && (bNext || (at[i].resind == resind)); i++)
122         {
123             if (anm[i] && gmx_strcasecmp(type, *(anm[i])) == 0)
124             {
125                 return i;
126             }
127         }
128         if (!(bNext && at[start].resind == at[natoms - 1].resind))
129         {
130             atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype,
131                            bAllowMissing);
132         }
133     }
134     else
135     {
136         /* The previous residue */
137         type++;
138         if (start > 0)
139         {
140             resind = at[start - 1].resind;
141         }
142         for (i = start - 1; (i >= 0) /*&& (at[i].resind == resind)*/; i--)
143         {
144             if (gmx_strcasecmp(type, *(anm[i])) == 0)
145             {
146                 return i;
147             }
148         }
149         if (start > 0)
150         {
151             atom_not_found(FARGS, type, at[start].resind, *atoms->resinfo[resind].name, bondtype,
152                            bAllowMissing);
153         }
154     }
155     return -1;
156 }
157
158 int search_res_atom(const char* type, int resind, const t_atoms* atoms, const char* bondtype, bool bAllowMissing)
159 {
160     int i;
161
162     for (i = 0; (i < atoms->nr); i++)
163     {
164         if (atoms->atom[i].resind == resind)
165         {
166             return search_atom(type, i, atoms, bondtype, bAllowMissing);
167         }
168     }
169
170     return -1;
171 }